corgidrp.astrom#
Functions#
|
Finds the center coordinates for a given frame. |
|
Finds the centroid in a sub-region around a given initial guess (or the brightest pixel if no guess is provided). |
|
Creates a template for matching psfs. |
|
Computes the relative offset between stars. |
|
Rough equivalivent to itertools.combinations function to create all r-length combinations from a given array. |
|
Used to find the angle from North counterclockwise between two sources in an image. |
|
Computes the linear distance between two points in polar coordinates. |
|
Converts position in separation (pixels from some reference center) and position angle |
|
Converts position in separation (pixels from some reference center) and position angle |
|
Used to find to [pixel, pixel] locations of the sources in an image |
|
Function to find the corresponding RA/Dec positions to image sources, given a particular field. |
|
Cost function used to fit the legendre polynomials for distortion mapping. |
|
Used to find the platescale and north angle of the image. Calculates the platescale for each pair of stars in the image |
|
Used to find the offset between the target and the center of the image. |
|
Function that formats the input data for the distortion map computation * must be run before compute_distortion * |
|
Function that computes the legendre polynomial coefficients that describe the image distortion map * must run format_disotrtio_inputs() first * |
|
Perform the boresight calibration of a dataset. |
|
Creates a circular mask |
|
Creates two, 2D maps of distortion for the X and Y directions from given legendre polynomial coefficients. |
Module Contents#
- corgidrp.astrom.centroid(frame)[source]#
Finds the center coordinates for a given frame.
- Parameters:
frame (np.ndarray) – 2D array to compute centering
- Returns:
xcen (float): X centroid coordinate ycen (float): Y centroid coordinate
- Return type:
tuple
- corgidrp.astrom.centroid_with_roi(frame, roi_radius=5, centering_initial_guess=None)[source]#
Finds the centroid in a sub-region around a given initial guess (or the brightest pixel if no guess is provided).
- Parameters:
frame (np.ndarray) – 2D array to compute centering.
roi_radius (int or float) – Half-size of the box around the initial guess or brightest pixel.
centering_initial_guess (tuple or None, optional) – (x_init, y_init) as initial guess for centroiding. If None, defaults to the brightest pixel.
- Returns:
xcen (float): X centroid coordinate. ycen (float): Y centroid coordinate.
- Return type:
tuple
- corgidrp.astrom.shift_psf(frame, dx, dy, flux, fitsize=10, stampsize=10)[source]#
Creates a template for matching psfs.
- Parameters:
frame (np.ndarray) – 2D array
dx (float) – Value to shift psf in ‘x’ direction
dy (float) – Value to shift psf in ‘y’ direction
flux (float) – Peak flux value
fitsize (float) – (Optional) Width of square on frame to fit psf to
stampsize (float) – (Optional) Width of psf stamp to create
- Returns:
New coordinates for a shifted psf
- Return type:
ndi..map_coordinates
- corgidrp.astrom.measure_offset(frame, xstar_guess, ystar_guess, xoffset_guess, yoffset_guess, guessflux=1, rad=5, stampsize=10)[source]#
Computes the relative offset between stars.
- Parameters:
frame (np.ndarray) – 2D array of data
xstar_guess (float) – Estimate of first star ‘x’ coordinate
ystar_guess (float) – Estimate of first star ‘y’ coordinate
xoffset_guess (float) – Estimate of ‘x’ direction offset between stars
yoffset_guess (float) – Estimate of ‘y’ direction offset between stars
guessflux (float) – (Optional) Peak flux of first star
rad (int) – (Optional) Radius around first star to compute centroid on
stampsize (float) – (Optional) Width of square psf stamp
- Returns:
List of [x,y] offsets in respective directions fit_errs (np.array): Array of [x,y] fitting errors
- Return type:
binary_offset (np.array)
- corgidrp.astrom.compute_combinations(iteration, r=2)[source]#
Rough equivalivent to itertools.combinations function to create all r-length combinations from a given array.
- Parameters:
iteration (np.array) – Array from which to create combinations of elements
r (int) – (Optional) Length of combinations (default: 2)
- corgidrp.astrom.angle_between(pos1, pos2)[source]#
Used to find the angle from North counterclockwise between two sources in an image.
- Parameters:
pos1 (tuple) – Position of the first target
pos2 (tuple) – Position of the second target
- Returns:
Angle [deg] between two sources, from north going counterclockwise
- Return type:
angle (float)
- corgidrp.astrom.get_polar_dist(seppa1, seppa2)[source]#
Computes the linear distance between two points in polar coordinates.
- Parameters:
seppa1 (tuple) – Separation (in any units) and position angle (in degrees) of the first point.
seppa2 (tuple) – Separation (in same units as above) and position angle (in degrees) of the second point.
- Returns:
Distance between the two points in the input separation units.
- Return type:
float
- corgidrp.astrom.seppa2dxdy(sep_pix, pa_deg)[source]#
Converts position in separation (pixels from some reference center) and position angle (counterclockwise from north) to separation in x and y pixels from the center.
- Parameters:
sep_pix (float or np.array) – Separation in pixels
pa_deg (float or np.array) – Position angle in degrees (counterclockwise from North)
- Returns:
array of shape (2,) containing delta x and delta y in pixels from the center
- Return type:
np.array
- corgidrp.astrom.seppa2xy(sep_pix, pa_deg, cenx, ceny)[source]#
Converts position in separation (pixels from some reference center) and position angle (counterclockwise from north) to separation in x and y pixels from the center.
- Parameters:
sep_pix (float or np.array) – Separation in pixels
pa_deg (float or np.array) – Position angle in degrees (counterclockwise from North)
cenx (float) – X location of center reference pixel. (0,0) is center of bottom left pixel
ceny (float) – Y location of center reference pixel. (0,0) is center of bottom left pixel
- Returns:
x and y pixel location. (0,0) is center of bottom left pixel.
- Return type:
np.array
- corgidrp.astrom.find_source_locations(image_data, threshold=10, fwhm=7, mask_rad=1)[source]#
Used to find to [pixel, pixel] locations of the sources in an image
- Parameters:
image_data (numpy.ndarray) – 2D array of image data
threshold (int) – Number of stars to find (default: 100)
fwhm (float) – Full width at half maximum of the stellar psf (default: 7, ~fwhm for a normal distribution with sigma=3)
mask_rad (int) – Radius of mask for stars [in fwhm] (default: 1)
- Returns:
Astropy table with columns ‘x’, ‘y’ as pixel locations
- Return type:
sources (astropy.table.Table)
- corgidrp.astrom.match_sources(image, sources, field_path, comparison_threshold=50, rad=0.012, platescale_guess=21.8, platescale_tol=0.1)[source]#
Function to find the corresponding RA/Dec positions to image sources, given a particular field.
- Parameters:
image (corgidrp.data.Image) – Image data as a corgidrp Image object
sources (astropy.table.Table) – Astropy table with columns ‘x’, ‘y’ as pixel locations of sources to match
field_path (str) – Full path to directory with search field data (ra, dec, vmag, etc.)
comparison_threshold (int) – How many stars in the field to consider for the initial match
rad (float) – The radius [deg] around the target coordinate for creating a subfield to match image sources to
platescale_guess (float) – An initial guess for the platescale value (default: 21.8 [mas/ pixel])
platescale_tol (float) – A tolerance for finding source matches within a fraction of the initial plate scale guess (default: 0.5)
- Returns:
Astropy table with columns ‘x’,’y’,’RA’, ‘DEC’ as pixel locations and corresponding sky positons
- Return type:
matched_sources (astropy.table.Table)
- corgidrp.astrom.fit_distortion_solution(params, fitorder, platescale, rotangle, pos1, meas_offset, sky_offset, meas_errs)[source]#
Cost function used to fit the legendre polynomials for distortion mapping.
- Parameters:
params (list) – List of the x and y legendre polynomial coefficients
fitorder (int) – The degree of legendre polynomial being used
platescale (float) – The platescale of the image
rotangle (float) – The north angle of the image
pos1 (np.array) – A (2 x N) array of (x, y) pixel positions for the first star in N pairs
meas_offset (np.array) – A (2 x N) array of (x, y) pixel offset from the first star position for N star pairs
sky_offset (np.array) – A (2 x N) array of (sep, pa) true sky offsets in [mas] and [deg] from the first star position for N pairs
meas_errs (np.array) – A (2 x N) array of (x, y) pixel errors in measured offsets from the first star position for N pairs
- Returns:
List of residuals between true and measured star positions
- Return type:
residuals (list)
- corgidrp.astrom.compute_platescale_and_northangle(image, source_info, center_coord, center_radius=0.9)[source]#
Used to find the platescale and north angle of the image. Calculates the platescale for each pair of stars in the image and returns the averged platescale. Calculates the north angle for pairs of stars with the center target and returns the averged north angle.
- Parameters:
image (numpy.ndarray) – 2D array of image data
source_info (astropy.table.Table) – Estimated pixel positions of sources and true sky positions, must have column names ‘x’, ‘y’, ‘RA’, ‘DEC’
center_coord (tuple) – (float): RA coordinate of the target pointing (float): Dec coordinate of the target pointing
center_radius (float) – Percent of the image radius used to crop the image and compute plate scale and north angle from (default: 1 – ie: the full image is used)
- Returns:
Platescale [mas/pixel] north_angle (float): Angle between image north and true north [deg]
- Return type:
platescale (float)
- corgidrp.astrom.compute_boresight(image, source_info, target_coordinate, cal_properties)[source]#
Used to find the offset between the target and the center of the image.
- Parameters:
image (numpy.ndarray) – 2D array of image data
source_info (astropy.table.Table) – Estimated pixel positions of sources and true sky positions, must have column names ‘x’, ‘y’, ‘RA’, ‘DEC’
target_coordinate (tuple) – (float): RA coordinate of the target pointing (float): DEC coordinate of the target pointing
cal_properties (tuple) – (float): Platescale (float): North angle
- Returns:
RA coordinate of the center pixel image_center_DEC (float): Dec coordinate of the center pixel
- Return type:
image_center_RA (float)
- corgidrp.astrom.format_distortion_inputs(input_dataset, source_matches, ref_star_pos, position_error=None)[source]#
Function that formats the input data for the distortion map computation * must be run before compute_distortion *
- Parameters:
input_dataset (corgidrp.data.dataset) – corgidrp dataset object with images to compute the distortion from
source_matches (list of astropy.table.Table() objects) – List of length N for N frames in the input dataset. Tables must columns ‘x’,’y’,’RA’,’DEC’ as pixel locations and corresponding sky positons
ref_star_pos (list of astropy.table.Table() objects) – List of length N for N frames. Tables must have column names ‘x’, ‘y’, ‘RA’, ‘DEC’ for the position of the reference position to compute pairs with
position_error (NoneType or int) – If int, this is the uniform error value assumed for the offset between pairs of stars in both x and y Should be changed later to accept non-uniform errors
- Returns:
2D array of the (x, y) pixel positions for the first star in every star pair offsets (np.array): 2D array of the (delta_x, delta_y) values for each star from the first star position true_offsets (np.array): 2D array of the (delta_ra, delta_dec) offsets between the matched stars in the reference field errs (np.array): 2D array of the (x_err, y_err) error in the measured pixel positions
- Return type:
first_stars (np.array)
- corgidrp.astrom.compute_distortion(input_dataset, pos1, meas_offset, sky_offset, meas_errs, platescale, northangle, fitorder=3, initial_guess=None)[source]#
Function that computes the legendre polynomial coefficients that describe the image distortion map * must run format_disotrtio_inputs() first *
- Parameters:
input_dataset (corgidrp.data.Dataset) – corgidrp dataset object with images to compute the distortion from
pos1 (np.array) – 2D array of the (x, y) pixel positions for the first star in every star pair
meas_offset (np.array) – 2D array of the (delta_x, delta_y) values for each star from the first star position
sky_offset (np.array) – 2D array of the (delta_ra, delta_dec) offsets between the matched stars in the reference field
meas_errs (np.array) – 2D array of the (x_err, y_err) error in the measured pixel positions
platescale (float) – Platescale value to use in computing distortion
northangle (float) – Northangle value to use in computing distortion
fitorder (int) – The order of legendre polynomial to fit to the image distortion (default: 3)
initial_guess (np.array) – Initial guess of fitting parameters (legendre coefficients) length based on fitorder (2 * (fitorder+1)**2), (default: None)
- Returns:
The legendre coefficients (np.array) and polynomial order used for the fit (int)
- Return type:
distortion_coeffs (tuple)
- corgidrp.astrom.boresight_calibration(input_dataset, field_path='JWST_CALFIELD2020.csv', field_matches=None, find_threshold=10, fwhm=7, mask_rad=1, comparison_threshold=50, search_rad=0.012, platescale_guess=21.8, platescale_tol=0.1, center_radius=0.9, frames_to_combine=None, find_distortion=False, fitorder=3, position_error=None, initial_dist_guess=None)[source]#
Perform the boresight calibration of a dataset.
- Parameters:
input_dataset (corgidrp.data.Dataset) – Dataset containing a images for astrometric calibration
field_path (str) – Full path to file with search field data (ra, dec, vmag, etc.) (default: ‘JWST_CALFIELD2020.csv’)
field_matches (list of str or astropy.table.Table) – List of full paths to files or astropy tables with calibration field matches for each image in the dataset (x, y, ra, dec), if single str the same filepath used for all frames,nif None, automated source matching is used (default: None)
find_threshold (int) – Number of stars to find (default 10)
fwhm (float) – Full width at half maximum of the stellar psf (default: 7, ~fwhm for a normal distribution with sigma=3)
mask_rad (int) – Radius of mask for stars [in fwhm] (default: 1)
comparison_threshold (int) – How many stars in the field to consider for the initial match (default: 50)
search_rad (float) – The radius [deg] around the target coordinate for creating a subfield to match image sources to
platescale_guess (float) – An initial guess for the platescale value (default: 21.8 [mas/ pixel])
platescale_tol (float) – A tolerance for finding source matches within a fraction of the initial plate scale guess (default: 0.1)
center_radius (float) – Percent of the image to compute plate scale and north angle from, centered around the image center (default: 0.9 – ie: 90% of the image is used)
frames_to_combine (int) – The number of frames to combine in a dataset (default: None)
find_distortion (boolean) – Used to determine if distortion map coeffs will be computed (default: False)
fitorder (int) – The order of legendre polynomials used to fit the distortion map (default: 3)
position_error (NoneType or int) – If int, this is the uniform error value assumed for the offset between pairs of stars in both x and y
initial_dist_guess (np.array) – An initial guess of legendre coefficients used for fitting distortion, if None will use coeffs associated with no distortion (default: None)
- Returns:
Astrometric Calibration data object containing image center coords in (RA,DEC), platescale, and north angle
- Return type:
- corgidrp.astrom.create_circular_mask(shape_yx, center=None, r=None)[source]#
Creates a circular mask
- Parameters:
shape_yx (list-like of int)
center (list of float, optional) – Center of mask. Defaults to the center of the array.
r (float, optional) – radius of mask. Defaults to the minimum distance from the center to the edge of the array.
- Returns:
boolean array with True inside the circle, False outside.
- Return type:
np.array
- corgidrp.astrom.transform_coeff_to_map(distortion_coeffs, fit_order, image_shape)[source]#
Creates two, 2D maps of distortion for the X and Y directions from given legendre polynomial coefficients.
- Parameters:
distortion_coeffs (array of float) – the array of legendre polynomial coefficients that describe the distortion map in x and y directions
fit_order (int) – the order of legendre polynomial used to fit distortion
image_shape (list-like of int) – the xy pixel shape of the image
- Returns:
a 2D array of distortion in the x direction across the image y_dist_map (np.ndarray): a 2D array of distortion in the y direction across the image
- Return type:
x_dist_map (np.ndarray)