Resample Utilities¶
- romancal.resample.resample_utils.build_driz_weight(model, weight_type=None, good_bits: str = None)[source]¶
Builds the drizzle weight map for resampling.
- Parameters:
model (object) – The input model.
weight_type (str, optional) – The type of weight to use. Allowed values are ‘ivm’ or ‘exptime’. Defaults to None.
good_bits (str, optional) – The good bits to use for building the mask. Defaults to None.
- Returns:
The drizzle weight map.
- Return type:
numpy.ndarray
- Raises:
ValueError – If an invalid weight type is provided.
Examples
model = get_input_model() weight_map = build_driz_weight(model, weight_type='ivm', good_bits=[1, 2, 3]) print(weight_map)
- romancal.resample.resample_utils.build_mask(dqarr, bitvalue)[source]¶
Build a bit mask from an input DQ array and a bitvalue flag.
- Parameters:
dqarr (numpy.ndarray) – Input DQ array.
bitvalue (str) – Bitvalue flag.
- Returns:
Bit mask where 1 represents good and 0 represents bad.
- Return type:
ndarray
Notes
The function interprets the bitvalue flag using the
astropy.nddata.bitmask.interpret_bit_flags
function.If the bitvalue is None, the function returns a bit mask with all elements set to 1.
Otherwise, the function performs a bitwise AND operation between the dqarr and the complement of the bitvalue, and then applies a logical NOT operation to obtain the bit mask.
The resulting bit mask is returned as a
numpy.ndarray
of dtypenumpy.uint8
.
- romancal.resample.resample_utils.calc_gwcs_pixmap(in_wcs, out_wcs, shape=None)[source]¶
Generate a pixel map grid using the input and output WCS.
- Parameters:
in_wcs (
WCS
) – Input WCS.out_wcs (
WCS
) – Output WCS.shape (tuple, optional) – Shape of the data. If provided, the bounding box will be calculated from the shape. If not provided, the bounding box will be calculated from the input WCS.
- Returns:
pixmap – The calculated pixel map grid.
- Return type:
ndarray
- romancal.resample.resample_utils.decode_context(context, x, y)[source]¶
Get 0-based indices of input images that contributed to (resampled) output pixel with coordinates
x
andy
.- Parameters:
context (numpy.ndarray) – A 3D
ndarray
of integral data type.x (int, list of integers, numpy.ndarray of integers) – X-coordinate of pixels to decode (3rd index into the
context
array)y (int, list of integers, numpy.ndarray of integers) – Y-coordinate of pixels to decode (2nd index into the
context
array)
- Returns:
A list of
numpy.ndarray
objects each containing indices of input imagesthat have contributed to an output pixel with coordinates
x
andy
.The length of returned list is equal to the number of input coordinate
arrays
x
andy
.
Examples
An example context array for an output image of array shape
(5, 6)
obtained by resampling 80 input images.con = np.array( [[[0, 0, 0, 0, 0, 0], [0, 0, 0, 36196864, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 537920000, 0, 0, 0]], [[0, 0, 0, 0, 0, 0,], [0, 0, 0, 67125536, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 163856, 0, 0, 0]], [[0, 0, 0, 0, 0, 0], [0, 0, 0, 8203, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 0, 0, 0, 0], [0, 0, 32865, 0, 0, 0]]], dtype=np.int32 ) decode_context(con, [3, 2], [1, 4]) [array([ 9, 12, 14, 19, 21, 25, 37, 40, 46, 58, 64, 65, 67, 77]), array([ 9, 20, 29, 36, 47, 49, 64, 69, 70, 79])]
- romancal.resample.resample_utils.make_output_wcs(input_models, pscale_ratio=None, pscale=None, rotation=None, shape=None, crpix: Tuple[float, float] = None, crval: Tuple[float, float] = None)[source]¶
Generate output WCS here based on footprints of all input WCS objects
- Parameters:
input_models (list of
roman_datamodels.datamodels.DataModel
) – Each datamodel must have agwcs.wcs.WCS
object.pscale_ratio (float, optional) – Ratio of input to output pixel scale. Ignored when
pscale
is provided.pscale (float, None, optional) – Absolute pixel scale in degrees. When provided, overrides
pscale_ratio
.rotation (float, None, optional) – Position angle (in degrees) of output image’s Y-axis relative to North. A value of 0.0 would orient the final output image to be North up. The default of
None
specifies that the images will not be rotated, but will instead be resampled in the default orientation for the camera with the x and y axes of the resampled image corresponding approximately to the detector axes.shape (tuple of int, None, optional) – Shape of the image (data array) using
numpy.ndarray
convention (ny
first andnx
second). This value will be assigned topixel_shape
andarray_shape
properties of the returned WCS object.crpix (tuple of float, None, optional) – Position of the reference pixel in the image array. If
ref_pixel
is not specified, it will be set to the center of the bounding box of the returned WCS object.crval (tuple of float, None, optional) – Right ascension and declination of the reference pixel. Automatically computed if not provided.
- Returns:
output_wcs – WCS object, with defined domain, covering entire set of input frames
- Return type:
object
- romancal.resample.resample_utils.reproject(wcs1, wcs2)[source]¶
Given two WCSs or transforms return a function which takes pixel coordinates in the first WCS or transform and computes them in the second one. It performs the forward transformation of
wcs1
followed by the inverse ofwcs2
.- Parameters:
wcs1 (
astropy.wcs.WCS
orgwcs.wcs.WCS
orastropy.modeling.Model
) – WCS objects.wcs2 (
astropy.wcs.WCS
orgwcs.wcs.WCS
orastropy.modeling.Model
) – WCS objects.
- Returns:
Function to compute the transformations. It takes
(x, y)
positions inwcs1
and returns(x, y)
positions inwcs2
.- Return type:
func