LSST Applications 24.1.5,g02d81e74bb+fa3a7a026e,g180d380827+a53a32eff8,g2079a07aa2+86d27d4dc4,g2305ad1205+c0501b3732,g295015adf3+7d3e92f0ec,g2bbee38e9b+0e5473021a,g337abbeb29+0e5473021a,g33d1c0ed96+0e5473021a,g3a166c0a6a+0e5473021a,g3ddfee87b4+5dd1654d75,g48712c4677+3bf1020dcb,g487adcacf7+065c13d9cf,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+d7ac436cfb,g5a732f18d5+53520f316c,g64a986408d+fa3a7a026e,g858d7b2824+fa3a7a026e,g8a8a8dda67+585e252eca,g99cad8db69+a5a909b84f,g9ddcbc5298+9a081db1e4,ga1e77700b3+15fc3df1f7,ga8c6da7877+4cf350ccb2,gb0e22166c9+60f28cb32d,gba4ed39666+c2a2e4ac27,gbb8dafda3b+f991a0b59f,gc120e1dc64+9ccbfdb8be,gc28159a63d+0e5473021a,gcf0d15dbbd+5dd1654d75,gd96a1ce819+42fd0ee607,gdaeeff99f8+f9a426f77a,ge6526c86ff+0d71447b4b,ge79ae78c31+0e5473021a,gee10cc3b42+585e252eca,gff1a9f87cc+fa3a7a026e
LSST Data Management Base Package
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB Class Reference

Public Member Functions

 __init__ (self, reference_array, log)
 
 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)
 

Public Attributes

 log
 

Protected Member Functions

 _build_distances_and_angles (self)
 
 _compute_test_vectors (self, source_array)
 
 _construct_pattern_and_shift_rot_matrix (self, src_pattern_array, n_match, max_cos_theta_shift, max_cos_rot_sq, max_dist_rad)
 
 _intermediate_verify (self, src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad)
 
 _create_spherical_rotation_matrix (self, rot_axis, cos_rotation, sin_rotion)
 
 _intermediate_verify_comparison (self, pattern_a, pattern_b, shift_rot_matrix, max_dist_rad)
 
 _test_pattern_lengths (self, test_pattern, max_dist_rad)
 
 _test_rotation_agreement (self, rot_vects, max_dist_rad)
 
 _final_verify (self, source_array, shift_rot_matrix, max_dist_rad, min_matches)
 
 _clip_distances (self, distances_rad)
 
 _match_sources (self, source_array, shift_rot_matrix)
 
 _handshake_match (self, matches_src, matches_ref)
 

Protected Attributes

 _reference_array
 
 _n_reference
 
 _dist_array
 
 _id_array
 

Detailed Description

Class implementing a pessimistic version of Optimistic Pattern Matcher
B (OPMb) from `Tabur (2007)`_, as described in `DMTN-031`_

Parameters
----------
reference_array : `numpy.ndarray`, (N, 3)
    spherical points x, y, z of to use as reference objects for
    pattern matching.
log : `lsst.log.Log` or `logging.Logger`
    Logger for outputting debug info.

Notes
-----
The class loads and stores the reference object
in a convenient data structure for matching any set of source objects that
are assumed to contain each other. The pessimistic nature of the algorithm
comes from requiring that it discovers at least two patterns that agree on
the correct shift and rotation for matching before exiting. The original
behavior of OPMb can be recovered simply. Patterns matched between the
input datasets are n-spoked pinwheels created from n+1 points. Refer to
`DMTN-031`_ for more details.

.. _Tabur (2007): https://doi.org/10.1071/AS07028
.. _DMTN-031: https://dmtn-031.lsst.io/

Definition at line 77 of file pessimistic_pattern_matcher_b_3D.py.

Constructor & Destructor Documentation

◆ __init__()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB.__init__ ( self,
reference_array,
log )

Definition at line 104 of file pessimistic_pattern_matcher_b_3D.py.

104 def __init__(self, reference_array, log):
105 self._reference_array = reference_array
106 self._n_reference = len(self._reference_array)
107 self.log = log
108
109 if self._n_reference > 0:
110 self._build_distances_and_angles()
111 else:
112 raise ValueError("No reference objects supplied")
113

Member Function Documentation

◆ _build_distances_and_angles()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._build_distances_and_angles ( self)
protected
Create the data structures we will use to search for our pattern
match in.

Throughout this function and the rest of the class we use id to
reference the position in the input reference catalog and index to
'index' into the arrays sorted on distance.

Definition at line 114 of file pessimistic_pattern_matcher_b_3D.py.

114 def _build_distances_and_angles(self):
115 """Create the data structures we will use to search for our pattern
116 match in.
117
118 Throughout this function and the rest of the class we use id to
119 reference the position in the input reference catalog and index to
120 'index' into the arrays sorted on distance.
121 """
122 # Create empty arrays to store our pair information per
123 # reference object.
124 self._dist_array = np.empty(
125 int(self._n_reference * (self._n_reference - 1) / 2),
126 dtype="float32")
127 self._id_array = np.empty(
128 (int(self._n_reference * (self._n_reference - 1) / 2), 2),
129 dtype="uint16")
130
131 startIdx = 0
132 # Loop over reference objects storing pair distances and ids.
133 for ref_id, ref_obj in enumerate(self._reference_array):
134 # Set the ending slicing index to the correct length for the
135 # pairs we are creating.
136 endIdx = startIdx + self._n_reference - 1 - ref_id
137
138 # Reserve and fill the ids of each reference object pair.
139 # 16 bit is safe for the id array as the catalog input from
140 # MatchPessimisticB is limited to a max length of 2 ** 16.
141 self._id_array[startIdx:endIdx, 0] = ref_id
142 self._id_array[startIdx:endIdx, 1] = np.arange(ref_id + 1,
143 self._n_reference,
144 dtype="uint16")
145
146 # Compute the vector deltas for each pair of reference objects.
147 # Compute and store the distances.
148 self._dist_array[startIdx:endIdx] = np.sqrt(
149 ((self._reference_array[ref_id + 1:, :]
150 - ref_obj) ** 2).sum(axis=1))
151 # Set startIdx of the slice to the end of the previous slice.
152 startIdx = endIdx
153
154 # Sort each array on the pair distances for the initial
155 # optimistic pattern matcher lookup.
156 sorted_dist_args = self._dist_array.argsort()
157 self._dist_array = self._dist_array[sorted_dist_args]
158 self._id_array = self._id_array[sorted_dist_args]
159

◆ _clip_distances()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._clip_distances ( self,
distances_rad )
protected
Compute a clipped max distance and calculate the number of pairs
that pass the clipped dist.

Parameters
----------
distances_rad : `numpy.ndarray`, (N,)
    Distances between pairs.

Returns
-------
output_struct : `lsst.pipe.base.Struct`
    Result struct with components:

    - ``n_matched_clipped`` : Number of pairs that survive the
      clipping on distance. (`float`)
    - ``clipped_max_dist`` : Maximum distance after clipping.
      (`float`).

Definition at line 764 of file pessimistic_pattern_matcher_b_3D.py.

764 def _clip_distances(self, distances_rad):
765 """Compute a clipped max distance and calculate the number of pairs
766 that pass the clipped dist.
767
768 Parameters
769 ----------
770 distances_rad : `numpy.ndarray`, (N,)
771 Distances between pairs.
772
773 Returns
774 -------
775 output_struct : `lsst.pipe.base.Struct`
776 Result struct with components:
777
778 - ``n_matched_clipped`` : Number of pairs that survive the
779 clipping on distance. (`float`)
780 - ``clipped_max_dist`` : Maximum distance after clipping.
781 (`float`).
782 """
783 clipped_dists, _, clipped_max_dist = sigmaclip(
784 distances_rad,
785 low=100,
786 high=2)
787 # Check clipped distances. The minimum value here
788 # prevents over convergence on perfect test data.
789 if clipped_max_dist < 1e-16:
790 clipped_max_dist = 1e-16
791 n_matched_clipped = np.sum(distances_rad < clipped_max_dist)
792 else:
793 n_matched_clipped = len(clipped_dists)
794
795 return pipeBase.Struct(n_matched_clipped=n_matched_clipped,
796 clipped_max_dist=clipped_max_dist)
797

◆ _compute_test_vectors()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._compute_test_vectors ( self,
source_array )
protected
Compute spherical 3 vectors at the edges of the x, y, z extent
of the input source catalog.

Parameters
----------
source_array : `numpy.ndarray`, (N, 3)
    array of 3 vectors representing positions on the unit
    sphere.

Returns
-------
test_vectors : `numpy.ndarray`, (N, 3)
    Array of vectors representing the maximum extents in x, y, z
    of the input source array. These are used with the rotations
    the code finds to test for agreement from different patterns
    when the code is running in pessimistic mode.

Definition at line 345 of file pessimistic_pattern_matcher_b_3D.py.

345 def _compute_test_vectors(self, source_array):
346 """Compute spherical 3 vectors at the edges of the x, y, z extent
347 of the input source catalog.
348
349 Parameters
350 ----------
351 source_array : `numpy.ndarray`, (N, 3)
352 array of 3 vectors representing positions on the unit
353 sphere.
354
355 Returns
356 -------
357 test_vectors : `numpy.ndarray`, (N, 3)
358 Array of vectors representing the maximum extents in x, y, z
359 of the input source array. These are used with the rotations
360 the code finds to test for agreement from different patterns
361 when the code is running in pessimistic mode.
362 """
363 # Get the center of source_array.
364 if np.any(np.logical_not(np.isfinite(source_array))):
365 self.log.warning("Input source objects contain non-finite values. "
366 "This could end badly.")
367 center_vect = np.nanmean(source_array, axis=0)
368
369 # So that our rotation test works over the full sky we compute
370 # the max extent in each Cartesian direction x,y,z.
371 xbtm_vect = np.array([np.min(source_array[:, 0]), center_vect[1],
372 center_vect[2]], dtype=np.float64)
373 xtop_vect = np.array([np.max(source_array[:, 0]), center_vect[1],
374 center_vect[2]], dtype=np.float64)
375 xbtm_vect /= np.sqrt(np.dot(xbtm_vect, xbtm_vect))
376 xtop_vect /= np.sqrt(np.dot(xtop_vect, xtop_vect))
377
378 ybtm_vect = np.array([center_vect[0], np.min(source_array[:, 1]),
379 center_vect[2]], dtype=np.float64)
380 ytop_vect = np.array([center_vect[0], np.max(source_array[:, 1]),
381 center_vect[2]], dtype=np.float64)
382 ybtm_vect /= np.sqrt(np.dot(ybtm_vect, ybtm_vect))
383 ytop_vect /= np.sqrt(np.dot(ytop_vect, ytop_vect))
384
385 zbtm_vect = np.array([center_vect[0], center_vect[1],
386 np.min(source_array[:, 2])], dtype=np.float64)
387 ztop_vect = np.array([center_vect[0], center_vect[1],
388 np.max(source_array[:, 2])], dtype=np.float64)
389 zbtm_vect /= np.sqrt(np.dot(zbtm_vect, zbtm_vect))
390 ztop_vect /= np.sqrt(np.dot(ztop_vect, ztop_vect))
391
392 # Return our list of vectors for later rotation testing.
393 return np.array([xbtm_vect, xtop_vect, ybtm_vect, ytop_vect,
394 zbtm_vect, ztop_vect])
395

◆ _construct_pattern_and_shift_rot_matrix()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._construct_pattern_and_shift_rot_matrix ( self,
src_pattern_array,
n_match,
max_cos_theta_shift,
max_cos_rot_sq,
max_dist_rad )
protected
Test an input source pattern against the reference catalog.
Returns the candidate matched patterns and their
implied rotation matrices or None.

Parameters
----------
src_pattern_array : `numpy.ndarray`, (N, 3)
    Sub selection of source 3 vectors to create a pattern from
n_match : `int`
    Number of points to attempt to create a pattern from. Must be
    >= len(src_pattern_array)
max_cos_theta_shift : `float`
    Maximum shift allowed between two patterns' centers.
max_cos_rot_sq : `float`
    Maximum rotation between two patterns that have been shifted
    to have their centers on top of each other.
max_dist_rad : `float`
    Maximum delta distance allowed between the source and reference
    pair distances to consider the reference pair a candidate for
    the source pair. Also sets the tolerance between the opening
    angles of the spokes when compared to the reference.

Return
-------
output_matched_pattern : `lsst.pipe.base.Struct`
    Result struct with components:

    - ``ref_candidates`` : ids of the matched pattern in the internal
      reference_array object (`list` of `int`).
    - ``src_candidates`` : Pattern ids of the sources matched
      (`list` of `int`).
    - ``shift_rot_matrix_src_to_ref`` : 3x3 matrix specifying the full
      shift and rotation between the reference and source objects.
      Rotates source into reference frame. `None` if match is not
      found. (`numpy.ndarray`, (3, 3))
    - ``shift_rot_matrix_ref_to_src`` : 3x3 matrix specifying the full
      shift and rotation of the reference and source objects. Rotates
      reference into source frame. None if match is not found
      (`numpy.ndarray`, (3, 3)).
    - ``cos_shift`` : Magnitude of the shift found between the two
      patten centers. `None` if match is not found (`float`).
    - ``sin_rot`` : float value of the rotation to align the already
      shifted source pattern to the reference pattern. `None` if no match
      found (`float`).

Definition at line 396 of file pessimistic_pattern_matcher_b_3D.py.

398 max_cos_rot_sq, max_dist_rad):
399 """Test an input source pattern against the reference catalog.
400 Returns the candidate matched patterns and their
401 implied rotation matrices or None.
402
403 Parameters
404 ----------
405 src_pattern_array : `numpy.ndarray`, (N, 3)
406 Sub selection of source 3 vectors to create a pattern from
407 n_match : `int`
408 Number of points to attempt to create a pattern from. Must be
409 >= len(src_pattern_array)
410 max_cos_theta_shift : `float`
411 Maximum shift allowed between two patterns' centers.
412 max_cos_rot_sq : `float`
413 Maximum rotation between two patterns that have been shifted
414 to have their centers on top of each other.
415 max_dist_rad : `float`
416 Maximum delta distance allowed between the source and reference
417 pair distances to consider the reference pair a candidate for
418 the source pair. Also sets the tolerance between the opening
419 angles of the spokes when compared to the reference.
420
421 Return
422 -------
423 output_matched_pattern : `lsst.pipe.base.Struct`
424 Result struct with components:
425
426 - ``ref_candidates`` : ids of the matched pattern in the internal
427 reference_array object (`list` of `int`).
428 - ``src_candidates`` : Pattern ids of the sources matched
429 (`list` of `int`).
430 - ``shift_rot_matrix_src_to_ref`` : 3x3 matrix specifying the full
431 shift and rotation between the reference and source objects.
432 Rotates source into reference frame. `None` if match is not
433 found. (`numpy.ndarray`, (3, 3))
434 - ``shift_rot_matrix_ref_to_src`` : 3x3 matrix specifying the full
435 shift and rotation of the reference and source objects. Rotates
436 reference into source frame. None if match is not found
437 (`numpy.ndarray`, (3, 3)).
438 - ``cos_shift`` : Magnitude of the shift found between the two
439 patten centers. `None` if match is not found (`float`).
440 - ``sin_rot`` : float value of the rotation to align the already
441 shifted source pattern to the reference pattern. `None` if no match
442 found (`float`).
443 """
444 # Create the delta vectors and distances we will need to assemble the
445 # spokes of the pattern.
446 src_delta_array = np.empty((len(src_pattern_array) - 1, 3))
447 src_delta_array[:, 0] = (src_pattern_array[1:, 0]
448 - src_pattern_array[0, 0])
449 src_delta_array[:, 1] = (src_pattern_array[1:, 1]
450 - src_pattern_array[0, 1])
451 src_delta_array[:, 2] = (src_pattern_array[1:, 2]
452 - src_pattern_array[0, 2])
453 src_dist_array = np.sqrt(src_delta_array[:, 0]**2
454 + src_delta_array[:, 1]**2
455 + src_delta_array[:, 2]**2)
456
457 pattern_result = construct_pattern_and_shift_rot_matrix(
458 src_pattern_array, src_delta_array, src_dist_array,
459 self._dist_array, self._id_array, self._reference_array, n_match,
460 max_cos_theta_shift, max_cos_rot_sq, max_dist_rad)
461
462 if pattern_result.success:
463 candidate_array = np.array(pattern_result.candidate_pairs)
464 fit_shift_rot_matrix = self._intermediate_verify(
465 src_pattern_array[candidate_array[:, 1]],
466 self._reference_array[candidate_array[:, 0]],
467 pattern_result.shift_rot_matrix, max_dist_rad)
468 return pipeBase.Struct(ref_candidates=candidate_array[:, 0].tolist(),
469 src_candidates=candidate_array[:, 1].tolist(),
470 shift_rot_matrix=fit_shift_rot_matrix,
471 cos_shift=pattern_result.cos_shift,
472 sin_rot=pattern_result.sin_rot)
473 return pipeBase.Struct(ref_candidates=[],
474 src_candidates=[],
475 shift_rot_matrix=None,
476 cos_shift=None,
477 sin_rot=None)
478

◆ _create_spherical_rotation_matrix()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._create_spherical_rotation_matrix ( self,
rot_axis,
cos_rotation,
sin_rotion )
protected
Construct a generalized 3D rotation matrix about a given
axis.

Parameters
----------
rot_axis : `numpy.ndarray`, (3,)
    3 vector defining the axis to rotate about.
cos_rotation : `float`
    cosine of the rotation angle.
sin_rotation : `float`
    sine of the rotation angle.

Return
------
shift_matrix : `numpy.ndarray`, (3, 3)
    3x3 spherical, rotation matrix.

Definition at line 535 of file pessimistic_pattern_matcher_b_3D.py.

536 sin_rotion):
537 """Construct a generalized 3D rotation matrix about a given
538 axis.
539
540 Parameters
541 ----------
542 rot_axis : `numpy.ndarray`, (3,)
543 3 vector defining the axis to rotate about.
544 cos_rotation : `float`
545 cosine of the rotation angle.
546 sin_rotation : `float`
547 sine of the rotation angle.
548
549 Return
550 ------
551 shift_matrix : `numpy.ndarray`, (3, 3)
552 3x3 spherical, rotation matrix.
553 """
554 rot_cross_matrix = np.array(
555 [[0., -rot_axis[2], rot_axis[1]],
556 [rot_axis[2], 0., -rot_axis[0]],
557 [-rot_axis[1], rot_axis[0], 0.]], dtype=np.float64)
558 shift_matrix = (cos_rotation*np.identity(3)
559 + sin_rotion*rot_cross_matrix
560 + (1. - cos_rotation)*np.outer(rot_axis, rot_axis))
561
562 return shift_matrix
563

◆ _final_verify()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._final_verify ( self,
source_array,
shift_rot_matrix,
max_dist_rad,
min_matches )
protected
Match the all sources into the reference catalog using the shift/rot
matrix.

After the initial shift/rot matrix is found, we refit the shift/rot
matrix using the matches the initial matrix produces to find a more
stable solution.

Parameters
----------
source_array : `numpy.ndarray` (N, 3)
    3-vector positions on the unit-sphere representing the sources to
    match
shift_rot_matrix : `numpy.ndarray` (3, 3)
    Rotation matrix representing inferred shift/rotation of the
    sources onto the reference catalog. Matrix need not be unitary.
max_dist_rad : `float`
    Maximum distance allowed for a match.
min_matches : `int`
    Minimum number of matched objects required to consider the
    match good.

Returns
-------
output_struct : `lsst.pipe.base.Struct`
    Result struct with components:

    - ``match_ids`` : Pairs of indexes into the source and reference
      data respectively defining a match (`numpy.ndarray`, (N, 2)).
    - ``distances_rad`` : distances to between the matched objects in
      the shift/rotated frame. (`numpy.ndarray`, (N,)).
    - ``max_dist_rad`` : Value of the max matched distance. Either
      returning the input value of the 2 sigma clipped value of the
      shift/rotated distances. (`float`).

Definition at line 665 of file pessimistic_pattern_matcher_b_3D.py.

669 min_matches):
670 """Match the all sources into the reference catalog using the shift/rot
671 matrix.
672
673 After the initial shift/rot matrix is found, we refit the shift/rot
674 matrix using the matches the initial matrix produces to find a more
675 stable solution.
676
677 Parameters
678 ----------
679 source_array : `numpy.ndarray` (N, 3)
680 3-vector positions on the unit-sphere representing the sources to
681 match
682 shift_rot_matrix : `numpy.ndarray` (3, 3)
683 Rotation matrix representing inferred shift/rotation of the
684 sources onto the reference catalog. Matrix need not be unitary.
685 max_dist_rad : `float`
686 Maximum distance allowed for a match.
687 min_matches : `int`
688 Minimum number of matched objects required to consider the
689 match good.
690
691 Returns
692 -------
693 output_struct : `lsst.pipe.base.Struct`
694 Result struct with components:
695
696 - ``match_ids`` : Pairs of indexes into the source and reference
697 data respectively defining a match (`numpy.ndarray`, (N, 2)).
698 - ``distances_rad`` : distances to between the matched objects in
699 the shift/rotated frame. (`numpy.ndarray`, (N,)).
700 - ``max_dist_rad`` : Value of the max matched distance. Either
701 returning the input value of the 2 sigma clipped value of the
702 shift/rotated distances. (`float`).
703 """
704 output_struct = pipeBase.Struct(
705 match_ids=None,
706 distances_rad=None,
707 max_dist_rad=None,
708 )
709
710 # Perform an iterative final verify.
711 match_sources_struct = self._match_sources(source_array,
712 shift_rot_matrix)
713 cut_ids = match_sources_struct.match_ids[
714 match_sources_struct.distances_rad < max_dist_rad]
715
716 n_matched = len(cut_ids)
717 clipped_struct = self._clip_distances(
718 match_sources_struct.distances_rad)
719 n_matched_clipped = clipped_struct.n_matched_clipped
720
721 if n_matched < min_matches or n_matched_clipped < min_matches:
722 return output_struct
723
724 # The shift/rotation matrix returned by
725 # ``_construct_pattern_and_shift_rot_matrix``, above, was
726 # based on only six points. Here, we refine that result by
727 # using all of the good matches from the “final verification”
728 # step, above. This will produce a more consistent result.
729 fit_shift_rot_matrix = least_squares(
730 _rotation_matrix_chi_sq,
731 x0=shift_rot_matrix.flatten(),
732 args=(source_array[cut_ids[:, 0], :3],
733 self._reference_array[cut_ids[:, 1], :3],
734 max_dist_rad)
735 ).x.reshape((3, 3))
736
737 # Redo the matching using the newly fit shift/rotation matrix.
738 match_sources_struct = self._match_sources(
739 source_array, fit_shift_rot_matrix)
740
741 # Double check the match distances to make sure enough matches
742 # survive still. We'll just overwrite the previously used variables.
743 n_matched = np.sum(
744 match_sources_struct.distances_rad < max_dist_rad)
745 clipped_struct = self._clip_distances(
746 match_sources_struct.distances_rad)
747 n_matched_clipped = clipped_struct.n_matched_clipped
748 clipped_max_dist = clipped_struct.clipped_max_dist
749
750 if n_matched < min_matches or n_matched_clipped < min_matches:
751 return output_struct
752
753 # Store our matches in the output struct. Decide between the clipped
754 # distance and the input max dist based on which is smaller.
755 output_struct.match_ids = match_sources_struct.match_ids
756 output_struct.distances_rad = match_sources_struct.distances_rad
757 if clipped_max_dist < max_dist_rad:
758 output_struct.max_dist_rad = clipped_max_dist
759 else:
760 output_struct.max_dist_rad = max_dist_rad
761
762 return output_struct
763

◆ _handshake_match()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._handshake_match ( self,
matches_src,
matches_ref )
protected
Return only those matches where both the source
and reference objects agree they they are each others'
nearest neighbor.

Parameters
----------
matches_src : `numpy.ndarray`, (N, 2)
    int array of nearest neighbor matches between shifted and
    rotated reference objects matched into the sources.
matches_ref : `numpy.ndarray`, (N, 2)
    int array of nearest neighbor matches between shifted and
    rotated source objects matched into the references.

Return
------
handshake_mask_array : `numpy.ndarray`, (N,)
    Array positions where the two match catalogs agree.

Definition at line 862 of file pessimistic_pattern_matcher_b_3D.py.

862 def _handshake_match(self, matches_src, matches_ref):
863 """Return only those matches where both the source
864 and reference objects agree they they are each others'
865 nearest neighbor.
866
867 Parameters
868 ----------
869 matches_src : `numpy.ndarray`, (N, 2)
870 int array of nearest neighbor matches between shifted and
871 rotated reference objects matched into the sources.
872 matches_ref : `numpy.ndarray`, (N, 2)
873 int array of nearest neighbor matches between shifted and
874 rotated source objects matched into the references.
875
876 Return
877 ------
878 handshake_mask_array : `numpy.ndarray`, (N,)
879 Array positions where the two match catalogs agree.
880 """
881 handshake_mask_array = np.zeros(len(matches_src), dtype=bool)
882
883 for src_match_idx, match in enumerate(matches_src):
884 ref_match_idx = np.searchsorted(matches_ref[:, 1], match[1])
885 if match[0] == matches_ref[ref_match_idx, 0]:
886 handshake_mask_array[src_match_idx] = True
887 return handshake_mask_array

◆ _intermediate_verify()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._intermediate_verify ( self,
src_pattern,
ref_pattern,
shift_rot_matrix,
max_dist_rad )
protected
 Perform an intermediate verify step.
Rotate the matches references into the source frame and test their
distances against tolerance. Only return true if all points are within
tolerance.

Parameters
----------
src_pattern : `numpy.ndarray`, (N,3)
    Array of 3 vectors representing the points that make up our source
    pinwheel pattern.
ref_pattern : `numpy.ndarray`, (N,3)
    Array of 3 vectors representing our candidate reference pinwheel
    pattern.
shift_rot_matrix : `numpy.ndarray`, (3,3)
    3x3 rotation matrix that takes the source objects and rotates them
    onto the frame of the reference objects
max_dist_rad : `float`
    Maximum distance allowed to consider two objects the same.

Returns
-------
fit_shift_rot_matrix : `numpy.ndarray`, (3,3)
   Fitted shift/rotation matrix if all of the points in our
   source pattern are within max_dist_rad of their matched reference
   objects. Returns None if this criteria is not satisfied.

Definition at line 479 of file pessimistic_pattern_matcher_b_3D.py.

480 max_dist_rad):
481 """ Perform an intermediate verify step.
482 Rotate the matches references into the source frame and test their
483 distances against tolerance. Only return true if all points are within
484 tolerance.
485
486 Parameters
487 ----------
488 src_pattern : `numpy.ndarray`, (N,3)
489 Array of 3 vectors representing the points that make up our source
490 pinwheel pattern.
491 ref_pattern : `numpy.ndarray`, (N,3)
492 Array of 3 vectors representing our candidate reference pinwheel
493 pattern.
494 shift_rot_matrix : `numpy.ndarray`, (3,3)
495 3x3 rotation matrix that takes the source objects and rotates them
496 onto the frame of the reference objects
497 max_dist_rad : `float`
498 Maximum distance allowed to consider two objects the same.
499
500 Returns
501 -------
502 fit_shift_rot_matrix : `numpy.ndarray`, (3,3)
503 Fitted shift/rotation matrix if all of the points in our
504 source pattern are within max_dist_rad of their matched reference
505 objects. Returns None if this criteria is not satisfied.
506 """
507 if len(src_pattern) != len(ref_pattern):
508 raise ValueError(
509 "Source pattern length does not match ref pattern.\n"
510 "\t source pattern len=%i, reference pattern len=%i" %
511 (len(src_pattern), len(ref_pattern)))
512
513 if self._intermediate_verify_comparison(
514 src_pattern, ref_pattern, shift_rot_matrix, max_dist_rad):
515 # Now that we know our initial shift and rot matrix is valid we
516 # want to fit the implied transform using all points from
517 # our pattern. This is a more robust rotation matrix as our
518 # initial matrix only used the first 2 points from the source
519 # pattern to estimate the shift and rotation. The matrix we fit
520 # are allowed to be non unitary but need to preserve the length of
521 # the rotated vectors to within the match tolerance.
522 fit_shift_rot_matrix = least_squares(
523 _rotation_matrix_chi_sq,
524 x0=shift_rot_matrix.flatten(),
525 args=(src_pattern, ref_pattern, max_dist_rad)
526 ).x.reshape((3, 3))
527 # Do another verify in case the fits have wondered off.
528 if self._intermediate_verify_comparison(
529 src_pattern, ref_pattern, fit_shift_rot_matrix,
530 max_dist_rad):
531 return fit_shift_rot_matrix
532
533 return None
534

◆ _intermediate_verify_comparison()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._intermediate_verify_comparison ( self,
pattern_a,
pattern_b,
shift_rot_matrix,
max_dist_rad )
protected
Test the input rotation matrix against one input pattern and
a second one.

If every point in the pattern after rotation is within a distance of
max_dist_rad to its candidate point in the other pattern, we return
True.

Parameters
----------
pattern_a :  `numpy.ndarray`, (N,3)
    Array of 3 vectors representing the points that make up our source
    pinwheel pattern.
pattern_b :  `numpy.ndarray`, (N,3)
    Array of 3 vectors representing our candidate reference pinwheel
    pattern.
shift_rot_matrix :  `numpy.ndarray`, (3,3)
    3x3 rotation matrix that takes the source objects and rotates them
    onto the frame of the reference objects
max_dist_rad : `float`
    Maximum distance allowed to consider two objects the same.


Returns
-------
bool
    True if all rotated source points are within max_dist_rad of
    the candidate references matches.

Definition at line 564 of file pessimistic_pattern_matcher_b_3D.py.

565 shift_rot_matrix, max_dist_rad):
566 """Test the input rotation matrix against one input pattern and
567 a second one.
568
569 If every point in the pattern after rotation is within a distance of
570 max_dist_rad to its candidate point in the other pattern, we return
571 True.
572
573 Parameters
574 ----------
575 pattern_a : `numpy.ndarray`, (N,3)
576 Array of 3 vectors representing the points that make up our source
577 pinwheel pattern.
578 pattern_b : `numpy.ndarray`, (N,3)
579 Array of 3 vectors representing our candidate reference pinwheel
580 pattern.
581 shift_rot_matrix : `numpy.ndarray`, (3,3)
582 3x3 rotation matrix that takes the source objects and rotates them
583 onto the frame of the reference objects
584 max_dist_rad : `float`
585 Maximum distance allowed to consider two objects the same.
586
587
588 Returns
589 -------
590 bool
591 True if all rotated source points are within max_dist_rad of
592 the candidate references matches.
593 """
594 shifted_pattern_a = np.dot(shift_rot_matrix,
595 pattern_a.transpose()).transpose()
596 tmp_delta_array = shifted_pattern_a - pattern_b
597 tmp_dist_array = (tmp_delta_array[:, 0] ** 2
598 + tmp_delta_array[:, 1] ** 2
599 + tmp_delta_array[:, 2] ** 2)
600 return np.all(tmp_dist_array < max_dist_rad ** 2)
601

◆ _match_sources()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._match_sources ( self,
source_array,
shift_rot_matrix )
protected
 Shift both the reference and source catalog to the the respective
frames and find their nearest neighbor using a kdTree.

Removes all matches who do not agree when either the reference or
source catalog is rotated. Cuts on a maximum distance are left to an
external function.

Parameters
----------
source_array : `numpy.ndarray`, (N, 3)
    array of 3 vectors representing the source objects we are trying
    to match into the source catalog.
shift_rot_matrix : `numpy.ndarray`, (3, 3)
    3x3 rotation matrix that performs the spherical rotation from the
    source frame into the reference frame.

Returns
-------
results : `lsst.pipe.base.Struct`
    Result struct with components:

    - ``matches`` : array of integer ids into the source and
      reference arrays. Matches are only returned for those that
      satisfy the distance and handshake criteria
      (`numpy.ndarray`, (N, 2)).
    - ``distances`` : Distances between each match in radians after
      the shift and rotation is applied (`numpy.ndarray`, (N)).

Definition at line 798 of file pessimistic_pattern_matcher_b_3D.py.

800 shift_rot_matrix):
801 """ Shift both the reference and source catalog to the the respective
802 frames and find their nearest neighbor using a kdTree.
803
804 Removes all matches who do not agree when either the reference or
805 source catalog is rotated. Cuts on a maximum distance are left to an
806 external function.
807
808 Parameters
809 ----------
810 source_array : `numpy.ndarray`, (N, 3)
811 array of 3 vectors representing the source objects we are trying
812 to match into the source catalog.
813 shift_rot_matrix : `numpy.ndarray`, (3, 3)
814 3x3 rotation matrix that performs the spherical rotation from the
815 source frame into the reference frame.
816
817 Returns
818 -------
819 results : `lsst.pipe.base.Struct`
820 Result struct with components:
821
822 - ``matches`` : array of integer ids into the source and
823 reference arrays. Matches are only returned for those that
824 satisfy the distance and handshake criteria
825 (`numpy.ndarray`, (N, 2)).
826 - ``distances`` : Distances between each match in radians after
827 the shift and rotation is applied (`numpy.ndarray`, (N)).
828 """
829 shifted_references = np.dot(
830 np.linalg.inv(shift_rot_matrix),
831 self._reference_array.transpose()).transpose()
832 shifted_sources = np.dot(
833 shift_rot_matrix,
834 source_array.transpose()).transpose()
835
836 ref_matches = np.empty((len(shifted_references), 2),
837 dtype="uint16")
838 src_matches = np.empty((len(shifted_sources), 2),
839 dtype="uint16")
840
841 ref_matches[:, 1] = np.arange(len(shifted_references),
842 dtype="uint16")
843 src_matches[:, 0] = np.arange(len(shifted_sources),
844 dtype="uint16")
845
846 ref_kdtree = cKDTree(self._reference_array)
847 src_kdtree = cKDTree(source_array)
848
849 ref_to_src_dist, tmp_ref_to_src_idx = \
850 src_kdtree.query(shifted_references)
851 src_to_ref_dist, tmp_src_to_ref_idx = \
852 ref_kdtree.query(shifted_sources)
853
854 ref_matches[:, 0] = tmp_ref_to_src_idx
855 src_matches[:, 1] = tmp_src_to_ref_idx
856
857 handshake_mask = self._handshake_match(src_matches, ref_matches)
858 return pipeBase.Struct(
859 match_ids=src_matches[handshake_mask],
860 distances_rad=src_to_ref_dist[handshake_mask],)
861

◆ _test_pattern_lengths()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._test_pattern_lengths ( self,
test_pattern,
max_dist_rad )
protected
 Test that the all vectors in a pattern are unit length within
tolerance.

This is useful for assuring the non unitary transforms do not contain
too much distortion.

Parameters
----------
test_pattern : `numpy.ndarray`, (N, 3)
    Test vectors at the maximum and minimum x, y, z extents.
max_dist_rad : `float`
    maximum distance in radians to consider two points "agreeing" on
    a rotation

Returns
-------
test : `bool`
    Length tests pass.

Definition at line 602 of file pessimistic_pattern_matcher_b_3D.py.

602 def _test_pattern_lengths(self, test_pattern, max_dist_rad):
603 """ Test that the all vectors in a pattern are unit length within
604 tolerance.
605
606 This is useful for assuring the non unitary transforms do not contain
607 too much distortion.
608
609 Parameters
610 ----------
611 test_pattern : `numpy.ndarray`, (N, 3)
612 Test vectors at the maximum and minimum x, y, z extents.
613 max_dist_rad : `float`
614 maximum distance in radians to consider two points "agreeing" on
615 a rotation
616
617 Returns
618 -------
619 test : `bool`
620 Length tests pass.
621 """
622 dists = (test_pattern[:, 0] ** 2
623 + test_pattern[:, 1] ** 2
624 + test_pattern[:, 2] ** 2)
625 return np.all(
626 np.logical_and((1 - max_dist_rad) ** 2 < dists,
627 dists < (1 + max_dist_rad) ** 2))
628

◆ _test_rotation_agreement()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._test_rotation_agreement ( self,
rot_vects,
max_dist_rad )
protected
 Test this rotation against the previous N found and return
the number that a agree within tolerance to where our test
points are.

Parameters
----------
rot_vects : `numpy.ndarray`, (N, 3)
    Arrays of rotated 3 vectors representing the maximum x, y,
    z extent on the unit sphere of the input source objects rotated by
    the candidate rotations into the reference frame.
max_dist_rad : `float`
    maximum distance in radians to consider two points "agreeing" on
    a rotation

Returns
-------
tot_consent : `int`
    Number of candidate rotations that agree for all of the rotated
    test 3 vectors.

Definition at line 629 of file pessimistic_pattern_matcher_b_3D.py.

629 def _test_rotation_agreement(self, rot_vects, max_dist_rad):
630 """ Test this rotation against the previous N found and return
631 the number that a agree within tolerance to where our test
632 points are.
633
634 Parameters
635 ----------
636 rot_vects : `numpy.ndarray`, (N, 3)
637 Arrays of rotated 3 vectors representing the maximum x, y,
638 z extent on the unit sphere of the input source objects rotated by
639 the candidate rotations into the reference frame.
640 max_dist_rad : `float`
641 maximum distance in radians to consider two points "agreeing" on
642 a rotation
643
644 Returns
645 -------
646 tot_consent : `int`
647 Number of candidate rotations that agree for all of the rotated
648 test 3 vectors.
649 """
650 self.log.debug("Comparing pattern %i to previous %i rotations...",
651 rot_vects[-1][-1], len(rot_vects) - 1)
652
653 tot_consent = 0
654 for rot_idx in range(max((len(rot_vects) - 1), 0)):
655 tmp_dist_list = []
656 for vect_idx in range(len(rot_vects[rot_idx]) - 1):
657 tmp_delta_vect = (rot_vects[rot_idx][vect_idx]
658 - rot_vects[-1][vect_idx])
659 tmp_dist_list.append(
660 np.dot(tmp_delta_vect, tmp_delta_vect))
661 if np.all(np.array(tmp_dist_list) < max_dist_rad ** 2):
662 tot_consent += 1
663 return tot_consent
664
int max

◆ match()

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB.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 )
Match a given source catalog into the loaded reference catalog.

Given array of points on the unit sphere and tolerances, we
attempt to match a pinwheel like pattern between these input sources
and the reference objects this class was created with. This pattern
informs of the shift and rotation needed to align the input source
objects into the frame of the references.

Parameters
----------
source_array : `numpy.ndarray`, (N, 3)
    An array of spherical x,y,z coordinates and a magnitude in units
    of objects having a lower value for sorting. The array should be
    of shape (N, 4).
n_check  : `int`
    Number of sources to create a pattern from. Not all objects may be
    checked if n_match criteria is before looping through all n_check
    objects.
n_match : `int`
    Number of objects to use in constructing a pattern to match.
n_agree : `int`
    Number of found patterns that must agree on their shift and
    rotation before exiting. Set this value to 1 to recover the
    expected behavior of Optimistic Pattern Matcher B.
max_n_patters : `int`
    Number of patterns to create from the input source objects to
    attempt to match into the reference objects.
max_shift : `float`
    Maximum allowed shift to match patterns in arcseconds.
max_rotation : `float`
    Maximum allowed rotation between patterns in degrees.
max_dist : `float`
    Maximum distance in arcseconds allowed between candidate spokes in
    the source and reference objects. Also sets that maximum distance
    in the intermediate verify, pattern shift/rotation agreement, and
    final verify steps.
pattern_skip_array : `int`
    Patterns we would like to skip. This could be due to the pattern
    being matched on a previous iteration that we now consider invalid.
    This assumes the ordering of the source objects is the same
    between different runs of the matcher which, assuming no object
    has been inserted or the magnitudes have changed, it should be.

Returns
-------
output_struct : `lsst.pipe.base.Struct`
    Result struct with components

    - ``matches`` : (N, 2) array of matched ids for pairs. Empty list if no
      match found (`numpy.ndarray`, (N, 2) or `list`)
    - ``distances_rad`` : Radian distances between the matched objects.
      Empty list if no match found (`numpy.ndarray`, (N,))
    - ``pattern_idx``: Index of matched pattern. None if no match found
      (`int`).
    - ``shift`` : Magnitude for the shift between the source and reference
      objects in arcseconds. None if no match found (`float`).

Definition at line 160 of file pessimistic_pattern_matcher_b_3D.py.

162 min_matches, pattern_skip_array=None):
163 """Match a given source catalog into the loaded reference catalog.
164
165 Given array of points on the unit sphere and tolerances, we
166 attempt to match a pinwheel like pattern between these input sources
167 and the reference objects this class was created with. This pattern
168 informs of the shift and rotation needed to align the input source
169 objects into the frame of the references.
170
171 Parameters
172 ----------
173 source_array : `numpy.ndarray`, (N, 3)
174 An array of spherical x,y,z coordinates and a magnitude in units
175 of objects having a lower value for sorting. The array should be
176 of shape (N, 4).
177 n_check : `int`
178 Number of sources to create a pattern from. Not all objects may be
179 checked if n_match criteria is before looping through all n_check
180 objects.
181 n_match : `int`
182 Number of objects to use in constructing a pattern to match.
183 n_agree : `int`
184 Number of found patterns that must agree on their shift and
185 rotation before exiting. Set this value to 1 to recover the
186 expected behavior of Optimistic Pattern Matcher B.
187 max_n_patters : `int`
188 Number of patterns to create from the input source objects to
189 attempt to match into the reference objects.
190 max_shift : `float`
191 Maximum allowed shift to match patterns in arcseconds.
192 max_rotation : `float`
193 Maximum allowed rotation between patterns in degrees.
194 max_dist : `float`
195 Maximum distance in arcseconds allowed between candidate spokes in
196 the source and reference objects. Also sets that maximum distance
197 in the intermediate verify, pattern shift/rotation agreement, and
198 final verify steps.
199 pattern_skip_array : `int`
200 Patterns we would like to skip. This could be due to the pattern
201 being matched on a previous iteration that we now consider invalid.
202 This assumes the ordering of the source objects is the same
203 between different runs of the matcher which, assuming no object
204 has been inserted or the magnitudes have changed, it should be.
205
206 Returns
207 -------
208 output_struct : `lsst.pipe.base.Struct`
209 Result struct with components
210
211 - ``matches`` : (N, 2) array of matched ids for pairs. Empty list if no
212 match found (`numpy.ndarray`, (N, 2) or `list`)
213 - ``distances_rad`` : Radian distances between the matched objects.
214 Empty list if no match found (`numpy.ndarray`, (N,))
215 - ``pattern_idx``: Index of matched pattern. None if no match found
216 (`int`).
217 - ``shift`` : Magnitude for the shift between the source and reference
218 objects in arcseconds. None if no match found (`float`).
219 """
220 # Given our input source_array we sort on magnitude.
221 sorted_source_array = source_array[source_array[:, -1].argsort(), :3]
222 n_source = len(sorted_source_array)
223
224 # Initialize output struct.
225 output_match_struct = pipeBase.Struct(
226 match_ids=[],
227 distances_rad=[],
228 pattern_idx=None,
229 shift=None,
230 max_dist_rad=None,)
231
232 if n_source <= 0:
233 self.log.warning("Source object array is empty. Unable to match. Exiting matcher.")
234 return None
235
236 # To test if the shifts and rotations we find agree with each other,
237 # we first create two test points situated at the top and bottom of
238 # where the z axis on the sphere bisects the source catalog.
239 test_vectors = self._compute_test_vectors(source_array[:, :3])
240
241 # We now create an empty list of our resultant rotated vectors to
242 # compare the different rotations we find.
243 rot_vect_list = []
244
245 # Convert the tolerances to values we will use in the code.
246 max_cos_shift = np.cos(np.radians(max_shift / 3600.))
247 max_cos_rot_sq = np.cos(np.radians(max_rotation)) ** 2
248 max_dist_rad = np.radians(max_dist / 3600.)
249
250 # Loop through the sources from brightest to faintest, grabbing a
251 # chunk of n_check each time.
252 for pattern_idx in range(np.min((max_n_patterns,
253 n_source - n_match))):
254
255 # If this pattern is one that we matched on the past but we
256 # now want to skip, we do so here.
257 if pattern_skip_array is not None and \
258 np.any(pattern_skip_array == pattern_idx):
259 self.log.debug(
260 "Skipping previously matched bad pattern %i...",
261 pattern_idx)
262 continue
263 # Grab the sources to attempt to create this pattern.
264 pattern = sorted_source_array[
265 pattern_idx: np.min((pattern_idx + n_check, n_source)), :3]
266
267 # Construct a pattern given the number of points defining the
268 # pattern complexity. This is the start of the primary tests to
269 # match our source pattern into the reference objects.
270 construct_return_struct = \
271 self._construct_pattern_and_shift_rot_matrix(
272 pattern, n_match, max_cos_shift, max_cos_rot_sq,
273 max_dist_rad)
274
275 # Our struct is None if we could not match the pattern.
276 if construct_return_struct.ref_candidates is None or \
277 construct_return_struct.shift_rot_matrix is None or \
278 construct_return_struct.cos_shift is None or \
279 construct_return_struct.sin_rot is None:
280 continue
281
282 # Grab the output data from the Struct object.
283 ref_candidates = construct_return_struct.ref_candidates
284 shift_rot_matrix = construct_return_struct.shift_rot_matrix
285 cos_shift = construct_return_struct.cos_shift
286 sin_rot = construct_return_struct.sin_rot
287
288 # If we didn't match enough candidates we continue to the next
289 # pattern.
290 if len(ref_candidates) < n_match:
291 continue
292
293 # Now that we know our pattern and shift/rotation are valid we
294 # store the the rotated versions of our test points for later
295 # use.
296 tmp_rot_vect_list = []
297 for test_vect in test_vectors:
298 tmp_rot_vect_list.append(np.dot(shift_rot_matrix, test_vect))
299 # Since our test point are in the source frame, we can test if
300 # their lengths are mostly preserved under the transform.
301 if not self._test_pattern_lengths(np.array(tmp_rot_vect_list),
302 max_dist_rad):
303 continue
304
305 tmp_rot_vect_list.append(pattern_idx)
306 rot_vect_list.append(tmp_rot_vect_list)
307
308 # Test if we have enough rotations, which agree, or if we
309 # are in optimistic mode.
310 if self._test_rotation_agreement(rot_vect_list, max_dist_rad) < \
311 n_agree - 1:
312 continue
313
314 # Run the final verify step.
315 match_struct = self._final_verify(source_array[:, :3],
316 shift_rot_matrix,
317 max_dist_rad,
318 min_matches)
319 if match_struct.match_ids is None or \
320 match_struct.distances_rad is None or \
321 match_struct.max_dist_rad is None:
322 continue
323
324 # Convert the observed shift to arcseconds
325 shift = np.degrees(np.arccos(cos_shift)) * 3600.
326 # Print information to the logger.
327 self.log.debug("Succeeded after %i patterns.", pattern_idx)
328 self.log.debug("\tShift %.4f arcsec", shift)
329 self.log.debug("\tRotation: %.4f deg",
330 np.degrees(np.arcsin(sin_rot)))
331
332 # Fill the struct and return.
333 output_match_struct.match_ids = \
334 match_struct.match_ids
335 output_match_struct.distances_rad = \
336 match_struct.distances_rad
337 output_match_struct.pattern_idx = pattern_idx
338 output_match_struct.shift = shift
339 output_match_struct.max_dist_rad = match_struct.max_dist_rad
340 return output_match_struct
341
342 self.log.debug("Failed after %i patterns.", pattern_idx + 1)
343 return output_match_struct
344

Member Data Documentation

◆ _dist_array

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._dist_array
protected

Definition at line 124 of file pessimistic_pattern_matcher_b_3D.py.

◆ _id_array

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._id_array
protected

Definition at line 127 of file pessimistic_pattern_matcher_b_3D.py.

◆ _n_reference

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._n_reference
protected

Definition at line 106 of file pessimistic_pattern_matcher_b_3D.py.

◆ _reference_array

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB._reference_array
protected

Definition at line 105 of file pessimistic_pattern_matcher_b_3D.py.

◆ log

lsst.meas.astrom.pessimistic_pattern_matcher_b_3D.PessimisticPatternMatcherB.log

Definition at line 107 of file pessimistic_pattern_matcher_b_3D.py.


The documentation for this class was generated from the following file: