matcher.ciratefi — RST-BC Invariant template match algorithm¶
The matcher.ciratefi module
New in version 0.1.0.
- autocnet.matcher.ciratefi.cifi(template, search_image, thresh=90, use_percentile=True, radii=[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11], scales=[0.5, 0.57, 0.66, 0.76, 0.87, 1.0], verbose=False)[source]¶
Circular sampling filter (Cifi) uses projections of the template and search images on a set of circular rings to detect the first grade candidate pixels and the points’ corresponding best fit scales for Ciratefi.
A set of scales is applied to the template and is radially sampled for each radii ‘r’ passed in. The template sample is equal to sum of the grayscale values divided by 2*pi*r.
Each pixel in the search image is similarly sampled. Every pixel then gets correlated with the template samples at all scales. The scales with the highest correlation are the considered the ‘best fit scales’.
- Parameters
template (np.array) – The input search template used to ‘query’ the destination image
search_image (np.array) – The image or sub-image to be searched
thresh (float) – The correlation thresh hold above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
use_percentile (bool) – If True (default), thresh is a percentile instead of a hard strength value
radii (np.array) – The list of radii to use for radial sampling
scales (list) – The list of scales to be applied to the template image, best if a geometric series
verbose (bool) – Set to True in order to output images and text describing the outputs. Can cause a serious decrease in performance. False by default.
- Returns
fg_candidate_pixels (ndarray) – array of pixels that passed the filter in tuples (y,x)
best_scales (ndarray) – parrallel array of best scales for the first grade candidate points
- autocnet.matcher.ciratefi.ciratefi(template, search_image, upsampling=1, cifi_thresh=95, rafi_thresh=95, tefi_thresh=100, use_percentile=True, alpha=0.19634954084936207, scales=[0.5, 0.57, 0.66, 0.76, 0.87, 1.0], radii=[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11], verbose=False)[source]¶
- Parameters
template (ndarray) – The input search template used to ‘query’ the destination image
search_image (ndarray) – The image or sub-image to be searched
upsampling (int) – upsample degree
cifi_thresh (float) – The correlation thresh hold for cifi above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
rafi_thresh (float) – The correlation thresh hold for rafi above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
tefi_thresh (float) – The correlation thresh hold for tefi above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
use_percentile (bool) – If True (default), thresh is a percentile instead of a hard strength value
alpha (float) – Must be greater than 0, if alpha is greater than 2*pi, it is reduced to it’s equivalent angle in [0,2*pi]. alpha list = np.arange(0, 2*pi, alpha)
radii (np.array) – The list of radii to use for radial sampling
scales (list) – The list of scales to be applied to the template image, best if a geometric series
verbose (bool) – Set to True in order to output images and text describing the outputs. Can cause a serious decrease in performance. False by default.
- Returns
results – array of pixel in (y, x)
- Return type
ndarray
- autocnet.matcher.ciratefi.circ_mask(shape, center, radius)[source]¶
Generates a circular mask
- Parameters
shape (tuple) – tuple decribing the desired mask shape in (y,x)
center (tuple) – tuple describing the desired center for the circle
radius (float) – radius of circlular mask
- Returns
mask – circular mask of bools
- Return type
ndarray
- autocnet.matcher.ciratefi.radial_line_mask(shape, center, radius, alpha=0.19460421, atol=0.01)[source]¶
Generates a linear mask from center at angle alpha.
- Parameters
shape (tuple) – tuple decribing the desired mask shape in (y,x)
center (tuple) – tuple describing the desired center for the circle
radius (float) – radius of the line mask
alpha (float) – angle for the line mask
atol (float) – absolute tolerance for alpha, the higher the tolerance, the wider the angle bandwidth
- Returns
mask – linear mask of bools
- Return type
ndarray
- autocnet.matcher.ciratefi.rafi(template, search_image, candidate_pixels, best_scales, thresh=95, use_percentile=True, alpha=0.39269908169872414, radii=[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11], verbose=False)[source]¶
The seconds filter in Ciratefi, the Radial Sampling Filter (Rafi), uses projections of the template image and the search image on a set of radial lines to upgrade the first grade the candidate pixels from cefi to seconds grade candidate pixels along with there corresponding best fit rotation.
The template image is radially sampled at angles 0-2*pi at steps alpha and with the best fit radius (largest sampling radius from radii list that fits in the template image)
Sampling for each line equals the sum of the greyscales divided by the best fit radius.
The search image is similarly sampled at each candidate pixel and is correlated with the radial samples on the template. The best fit angle is the angle that maximizes this correlation, and the second grade candidate pixels are determined by the strength of the correlation and the passed threshold
- Parameters
template (np.array) – The input search template used to ‘query’ the destination image
search_image (np.array) – The image or sub-image to be searched
candidate_pixels (np.array) – array of candidate pixels in tuples (y,x), best if the pixel are the output of Cifi
best_scales (ndarray) – The list of best fit scales for each candidate point, the length should equal the length of the candidate point list
thresh (float) – The correlation thresh hold above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
use_percentile (bool) – If True (default), thresh is a percentile instead of a hard strength value
alpha (float) – Must be greater than 0, if alpha is greater than 2*pi, it is reduced to it’s equivalent angle in [0,2*pi]. alpha list = np.arange(0, 2*pi, alpha)
radii (list) – The list of radii to use for radial sampling, best if the list is the same as the one used for Cifi
verbose (bool) – Set to True in order to output images and text describing the outputs. Can cause a serious decrease in performance. False by default.
- Returns
sg_candidate_points (ndarray)
best_rotation (ndarray) – Parrallel array of the best fit rotations for each second grade candidate pixel
- autocnet.matcher.ciratefi.tefi(template, search_image, candidate_pixels, best_scales, best_angles, scales=[0.5, 0.57, 0.66, 0.76, 0.87, 1.0], upsampling=1, thresh=100, alpha=0.19634954084936207, use_percentile=True, verbose=False)[source]¶
Template Matching Filter (Tefi) is the third and final filter for ciratefi.
For every candidate pixel, tefi rotates and scales the template image by the list of scales and angles passed in (which, ideally are the output from cefi and rafi respectively) and performs template match around the candidate pixels at the approriate scale and rotation angle. Here, the scales, angles and candidate points should be a parrallel array structure.
Any points with correlation strength over the threshold are returned as the the strongest candidates for the image location. If knows the point exists in one location, thresh should be 100 and use_percentile = True.
- Parameters
template (ndarray) – The input search template used to ‘query’ the destination image
image (ndarray) – The image or sub-image to be searched
candidate_pixels (ndarray) – array of candidate pixels in tuples (y,x), best if the pixel are the output of Cifi
best_scales (ndarray) – The list of best fit scales for each candidate point, the length should equal the length of the candidate point list
best_angles (ndarray) – The list of best fit rotation for each candidate point in radians, the length should equal the length of the candidate point list
upsampling (int) – upsample degree
thresh (float) – The correlation thresh hold above which a point will be a first grade candidate point. If use_percentile=True this will act as a percentile, for example, passing 90 means keep values in the top 90th percentile
use_percentile (bool) – If True (default), thresh is a percentile instead of a hard strength value
alpha (float) – A float between 0 & 2*pi, alpha list = np.arange(0, 2*pi, alpha)
verbose (bool) – Set to True in order to output images and text describing the outputs. Can cause a serious decrease in performance. False by default.
- Returns
results – array of pixel tuples (y,x) which over the threshold
- Return type
ndarray
- autocnet.matcher.ciratefi.to_polar_coord(shape, center)[source]¶
Generate a polar coordinate grid from a shape given a center.
- Parameters
shape (tuple) – tuple decribing the desired shape in (y,x)
center (tuple) – tuple describing the desired center for the grid
- Returns
r2 (ndarray) – grid of radii from the center
theta (ndarray) – grid of angles from the center