1
1
# Licensed under a 3-clause BSD style license - see LICENSE.rst
2
2
3
3
"""This module implements cutout functionality similar to fitscut, but for the ASDF file format."""
4
+ import copy
4
5
import pathlib
5
- from typing import Union
6
+ from typing import Union , Tuple
6
7
7
8
import asdf
8
9
import astropy
9
10
import gwcs
10
11
import numpy as np
11
12
12
13
from astropy .coordinates import SkyCoord
14
+ from astropy .modeling import models
13
15
14
16
15
17
def get_center_pixel (gwcsobj : gwcs .wcs .WCS , ra : float , dec : float ) -> tuple :
@@ -56,7 +58,8 @@ def get_center_pixel(gwcsobj: gwcs.wcs.WCS, ra: float, dec: float) -> tuple:
56
58
57
59
def get_cutout (data : asdf .tags .core .ndarray .NDArrayType , coords : Union [tuple , SkyCoord ],
58
60
wcs : astropy .wcs .wcs .WCS = None , size : int = 20 , outfile : str = "example_roman_cutout.fits" ,
59
- write_file : bool = True , fill_value : Union [int , float ] = np .nan ) -> astropy .nddata .Cutout2D :
61
+ write_file : bool = True , fill_value : Union [int , float ] = np .nan ,
62
+ gwcsobj : gwcs .wcs .WCS = None ) -> astropy .nddata .Cutout2D :
60
63
""" Get a Roman image cutout
61
64
62
65
Cut out a square section from the input image data array. The ``coords`` can either be a tuple of x, y
@@ -79,6 +82,8 @@ def get_cutout(data: asdf.tags.core.ndarray.NDArrayType, coords: Union[tuple, Sk
79
82
Flag to write the cutout to a file or not
80
83
fill_value: int | float, by default np.nan
81
84
The fill value for pixels outside the original image.
85
+ gwcsobj : gwcs.wcs.WCS, Optional
86
+ the original gwcs object for the full image, needed only when writing cutout as asdf file
82
87
83
88
Returns
84
89
-------
@@ -91,6 +96,8 @@ def get_cutout(data: asdf.tags.core.ndarray.NDArrayType, coords: Union[tuple, Sk
91
96
when a wcs is not present when coords is a SkyCoord object
92
97
RuntimeError:
93
98
when the requested cutout does not overlap with the original image
99
+ ValueError:
100
+ when no gwcs object is provided when writing to an asdf file
94
101
"""
95
102
96
103
# check for correct inputs
@@ -122,12 +129,23 @@ def get_cutout(data: asdf.tags.core.ndarray.NDArrayType, coords: Union[tuple, Sk
122
129
if write_as == '.fits' :
123
130
_write_fits (cutout , outfile )
124
131
elif write_as == '.asdf' :
125
- _write_asdf (cutout , outfile )
132
+ if not gwcsobj :
133
+ raise ValueError ('The original gwcs object is needed when writing to asdf file.' )
134
+ _write_asdf (cutout , gwcsobj , outfile )
126
135
127
136
return cutout
128
137
129
138
130
- def _write_fits (cutout , outfile = "example_roman_cutout.fits" ):
139
+ def _write_fits (cutout : astropy .nddata .Cutout2D , outfile : str = "example_roman_cutout.fits" ):
140
+ """ Write cutout as FITS file
141
+
142
+ Parameters
143
+ ----------
144
+ cutout : astropy.nddata.Cutout2D
145
+ the 2d cutout
146
+ outfile : str, optional
147
+ the name of the output cutout file, by default "example_roman_cutout.fits"
148
+ """
131
149
# check if the data is a quantity and get the array data
132
150
if isinstance (cutout .data , astropy .units .Quantity ):
133
151
data = cutout .data .value
@@ -137,8 +155,61 @@ def _write_fits(cutout, outfile="example_roman_cutout.fits"):
137
155
astropy .io .fits .writeto (outfile , data = data , header = cutout .wcs .to_header (relax = True ), overwrite = True )
138
156
139
157
140
- def _write_asdf (cutout , outfile = "example_roman_cutout.asdf" ):
141
- tree = {'roman' : {'meta' : {'wcs' : dict (cutout .wcs .to_header (relax = True ))}, 'data' : cutout .data }}
158
+ def _slice_gwcs (gwcsobj : gwcs .wcs .WCS , slices : Tuple [slice , slice ]) -> gwcs .wcs .WCS :
159
+ """ Slice the original gwcs object
160
+
161
+ "Slices" the original gwcs object down to the cutout shape. This is a hack
162
+ until proper gwcs slicing is in place a la fits WCS slicing. The ``slices``
163
+ keyword input is a tuple with the x, y cutout boundaries in the original image
164
+ array, e.g. ``cutout.slices_original``. Astropy Cutout2D slices are in the form
165
+ ((ymin, ymax, None), (xmin, xmax, None))
166
+
167
+ Parameters
168
+ ----------
169
+ gwcsobj : gwcs.wcs.WCS
170
+ the original gwcs from the input image
171
+ slices : Tuple[slice, slice]
172
+ the cutout x, y slices as ((ymin, ymax), (xmin, xmax))
173
+
174
+ Returns
175
+ -------
176
+ gwcs.wcs.WCS
177
+ The sliced gwcs object
178
+ """
179
+ tmp = copy .deepcopy (gwcsobj )
180
+
181
+ # get the cutout array bounds and create a new shift transform to the cutout
182
+ # add the new transform to the gwcs
183
+ xmin , xmax = slices [1 ].start , slices [1 ].stop
184
+ ymin , ymax = slices [0 ].start , slices [0 ].stop
185
+ shape = (ymax - ymin , xmax - xmin )
186
+ offsets = models .Shift (xmin , name = 'cutout_offset1' ) & models .Shift (ymin , name = 'cutout_offset2' )
187
+ tmp .insert_transform ('detector' , offsets , after = True )
188
+
189
+ # modify the gwcs bounding box to the cutout shape
190
+ tmp .bounding_box = ((0 , shape [0 ] - 1 ), (0 , shape [1 ] - 1 ))
191
+ tmp .pixel_shape = shape [::- 1 ]
192
+ tmp .array_shape = shape
193
+ return tmp
194
+
195
+
196
+ def _write_asdf (cutout : astropy .nddata .Cutout2D , gwcsobj : gwcs .wcs .WCS , outfile : str = "example_roman_cutout.asdf" ):
197
+ """ Write cutout as ASDF file
198
+
199
+ Parameters
200
+ ----------
201
+ cutout : astropy.nddata.Cutout2D
202
+ the 2d cutout
203
+ gwcsobj : gwcs.wcs.WCS
204
+ the original gwcs object for the full image
205
+ outfile : str, optional
206
+ the name of the output cutout file, by default "example_roman_cutout.asdf"
207
+ """
208
+ # slice the origial gwcs to the cutout
209
+ sliced_gwcs = _slice_gwcs (gwcsobj , cutout .slices_original )
210
+
211
+ # create the asdf tree
212
+ tree = {'roman' : {'meta' : {'wcs' : sliced_gwcs , 'data' : cutout .data }}}
142
213
af = asdf .AsdfFile (tree )
143
214
144
215
# Write the data to a new file
@@ -186,4 +257,4 @@ def asdf_cut(input_file: str, ra: float, dec: float, cutout_size: int = 20,
186
257
187
258
# create the 2d image cutout
188
259
return get_cutout (data , pixel_coordinates , wcs , size = cutout_size , outfile = output_file ,
189
- write_file = write_file , fill_value = fill_value )
260
+ write_file = write_file , fill_value = fill_value , gwcsobj = gwcsobj )
0 commit comments