pypeit.core.wavecal.patterns module
Module for finding patterns in arc line spectra
- pypeit.core.wavecal.patterns.curved_quadrangles(detlines, linelist, npixels, detsrch=5, lstsrch=10, pixtol=1.0)[source]
Brute force pattern recognition using curved quadrangles. A curved quadrangle contains (for either detlines or linelist):
a left line (l),
a right line (r),
a mid line (m), and
one line in between (c; c != m)
Something like this:
> | > | | | > | | | | > | | | | > l c m r
Then, the values (c-r)/(r-l) are in the same coordinate system for both detlines and linelist.
- Parameters:
detlines (ndarray) – list of detected lines in pixels (sorted, increasing)
linelist (ndarray) – list of lines that should be detected (sorted, increasing)
npixels (float) – Number of pixels along the dispersion direction
detsrch (int) – Number of consecutive elements in detlines to use to create a pattern (-1 means all lines in detlines)
lstsrch (int) – Number of consecutive elements in linelist to use to create a pattern (-1 means all lines in detlines)
pixtol (float) – tolerance that is used to determine if a match is successful (in units of pixels)
- Returns:
dindex (ndarray) – Index array of all detlines used in each triangle
lindex (ndarray) – Index array of the assigned line to each index in dindex
wvcen (ndarray) – central wavelength of each triangle
disps (ndarray) – Dispersion of each triangle (angstroms/pixel)
- pypeit.core.wavecal.patterns.detect_2Dpeaks(image)[source]
Takes a 2D image and returns 1 if a local maximum is found, and 0 otherwise
- Parameters:
image (ndarray) – 2D image
- Returns:
pimage – boolean mask of the peaks (1 when the pixel’s value is the neighborhood maximum, 0 otherwise)
- Return type:
ndarray
- pypeit.core.wavecal.patterns.empty_patt_dict(nlines)[source]
Return an empty patt_dict
- Parameters:
nlines – Number of lines for creating the mask.
- Returns:
patt_dict – An empty pattern dictionary
- Return type:
- pypeit.core.wavecal.patterns.match_quad_to_list(spec_lines, line_list, wv_guess, dwv_guess, tol=2.0, dwv_uncertainty=0.2, min_ftol=0.005)[source]
- pypeit.core.wavecal.patterns.quadrangles(detlines, linelist, npixels, detsrch=5, lstsrch=10, pixtol=1.0)[source]
Brute force pattern recognition using quadrangles. A quadrangle contains (for either detlines or linelist):
a left line (l),
a right line (r), and
two lines in between (a, b)
Something like this:
> | > | | | > | | | | > | | | | > l a b r
Then, the values (a-ll)/(r-ll) and (b-ll)/(r-ll) are in the same coordinate system for both detlines and linelist.
- Parameters:
detlines (ndarray) – list of detected lines in pixels (sorted, increasing)
linelist (ndarray) – list of lines that should be detected (sorted, increasing)
npixels (float) – Number of pixels along the dispersion direction
detsrch (int) – Number of consecutive elements in detlines to use to create a pattern (-1 means all lines in detlines)
lstsrch (int) – Number of consecutive elements in linelist to use to create a pattern (-1 means all lines in detlines)
pixtol (float) – tolerance that is used to determine if a match is successful (in units of pixels)
- Returns:
dindex (ndarray) – Index array of all detlines used in each triangle
lindex (ndarray) – Index array of the assigned line to each index in dindex
wvcen (ndarray) – central wavelength of each triangle
disps (ndarray) – Dispersion of each triangle (angstroms/pixel)
- pypeit.core.wavecal.patterns.run_quad_match(tcent, twave, llist_wv, disp, swv_uncertainty=250.0, pix_tol=1.0)[source]
- Parameters:
tcent (ndarray) – Pixel positions of arc lines
twave (ndarray) – Crude guess at wavelength solution, e.g. from wvcen, disp
llist_wv (ndarray) – Lines to match against (from a line list)
pix_tol (float) – Tolerance in units of pixels to match to
- Returns:
match_idx (dict) – Record of matches
scores (ndarray) – str array of scores
- pypeit.core.wavecal.patterns.scan_for_matches(wvcen, disp, npix, cut_tcent, wvdata, best_dict=None, swv_uncertainty=350.0, wvoff=1000.0, pix_tol=2.0, ampl=None)[source]
Warning
best_dict is updated in place
- pypeit.core.wavecal.patterns.score_quad_matches(fidx)[source]
Grades quad_match results
- Parameters:
fidx –
- Returns:
scores
- Return type:
- pypeit.core.wavecal.patterns.score_triangles(counts)[source]
Grades for the triangle results
- Parameters:
counts (ndarray) – Each element is a counter, representing a wavelength that is attributed to a given detected line. The more times that a wavelength is attributed to a detected line, the higher the counts. The more different wavelengths that are attributed to the same detected line (i.e. not ideal) the longer the counts list will be.
- Returns:
score – A string indicating the relative quality of the ID
- Return type:
- pypeit.core.wavecal.patterns.score_xcorr(counts, cc_avg, nreid_min=4, cc_local_thresh=-1.0)[source]
Grades for the cross-correlation results
- Parameters:
counts (ndarray) – Each element is a counter, representing a wavelength that is attributed to a given detected line. The more times that a wavelength is attributed to a detected line, the higher the counts. The more different wavelengths that are attributed to the same detected line (i.e. not ideal) the longer the counts list will be.
nreid_min (int, default = 4, optional) – Minimum number of matches to receive a score of ‘Perfect’ or ‘Very Good’
cc_local_thresh (float, default = -1.0, optional) – What does this do??
- Returns:
score – A string indicating the relative quality of the ID
- Return type:
- pypeit.core.wavecal.patterns.solve_triangles(detlines, linelist, dindex, lindex, patt_dict=None)[source]
Given a starting solution, find the best match for all detlines
- Parameters:
detlines (ndarray) – list of detected lines in pixels (sorted, increasing)
linelist (ndarray) – list of lines that should be detected (sorted, increasing)
dindex (ndarray) – Index array of all detlines (pixels) used in each triangle
lindex (ndarray) – Index array of the assigned line (wavelengths)to each index in dindex
patt_dict (dict) – Contains all relevant details of the fit
- pypeit.core.wavecal.patterns.solve_xcorr(detlines, linelist, dindex, lindex, line_cc, nreid_min: int = 4, cc_local_thresh: float = 0.8)[source]
Given a starting solution, find the best match for all detlines
- Parameters:
detlines (numpy.ndarray) – list of detected lines in pixels (sorted, increasing)
linelist (numpy.ndarray) – list of lines that should be detected (sorted, increasing)
dindex (numpy.ndarray) – Index array of all detlines (pixels) used in each triangle
lindex (numpy.ndarray) – Index array of the assigned line (wavelengths)to each index in dindex
line_cc (numpy.ndarray) – local cross correlation coefficient computed for each line
cc_local_thresh (float, default = 0.8, optional) – Threshold to satisy for local cross-correlation
nreid_min (int, default = 4, optional) – Minimum number of matches to receive a score of ‘Perfect’ or ‘Very Good’ Passed to score_xcorr()
- Returns:
patt_dict –
Contains all relevant details of the IDs. Keys are:
acceptable: bool: flag indicating success or failure
mask: ndarray, dtype =bool: mask indicating which lines are good
nmatch: int: Number of matching lines
scores: ndarray, str: Scores of the lines
IDs: ndarray, float: Wavelength IDs of the lines
cc_avg: ndarray, float: Average local zero-lag cross-correlation (over all the spectra for which a match was obtained) for the most often occuring wavlength ID
- Return type:
- pypeit.core.wavecal.patterns.triangles(detlines, linelist, npixels, detsrch=5, lstsrch=10, pixtol=1.0)[source]
Brute force pattern recognition using triangles. A triangle contains (for either detlines or linelist):
a starting point (s),
an end point (e), and
something in between (b)
Something like this:
> | > | | > | | | > | | | > s b e
Then, the value (b-s)/(e-s) is in the same coordinate system for both detlines and linelist.
- Parameters:
detlines (ndarray) – list of detected lines in pixels (sorted, increasing)
linelist (ndarray) – list of lines that should be detected (sorted, increasing)
npixels (float) – Number of pixels along the dispersion direction
detsrch (int) – Number of consecutive elements in detlines to use to create a pattern (-1 means all lines in detlines)
lstsrch (int) – Number of consecutive elements in linelist to use to create a pattern (-1 means all lines in detlines)
pixtol (float) – tolerance that is used to determine if a match is successful (in units of pixels)
- Returns:
dindex (ndarray) – Index array of all detlines used in each triangle
lindex (ndarray) – Index array of the assigned line to each index in dindex
wvcen (ndarray) – central wavelength of each triangle
disps (ndarray) – Dispersion of each triangle (angstroms/pixel)