Skip to content

Commit

Permalink
update for cos(dec)
Browse files Browse the repository at this point in the history
  • Loading branch information
jemorrison committed Aug 29, 2024
1 parent b0011f5 commit 1c5a745
Show file tree
Hide file tree
Showing 2 changed files with 35 additions and 10 deletions.
19 changes: 10 additions & 9 deletions docs/jwst/cube_build/arguments.rst
Original file line number Diff line number Diff line change
Expand Up @@ -131,24 +131,25 @@ A parameter only used for investigating which detector pixels contributed to a c
To print information to the screeen about the x = 10, y = 20, z = 35 spaxel the parameter string value is '10 20 35'.

.. _offsets:
The offset file is an ASDF formated file : <https://asdf-standard.readthedocs.io/>`_ stands for "Advanced Scientific Data. For each
input file in the spec3 assocation used to build the IFU cubes, an ra and dec offset, in arc seconds, is provided.
The offset file is an ASDF formated file :`<https://asdf-standard.readthedocs.io/>`_ stands for "Advanced Scientific Data. For each
input file in the spec3 assocation used to build the IFU cubes, there is a corresponding right ascension and declination offset,
given arc seconds.
Below is an example of how to make an ASDF offset file. It is assumed the user has determined the
offsets to apply for each file. The offsets are stored in a python dictionary, `offsets`. The items of this dictionary are `filenames`, `raoffset` and `decoffset`. The cube_building code is expects this dictionary to hold the information
for storing the file names and the associated ra and dec offsets.
offsets to apply the data in each file. The offsets are stored in a python dictionary, `offsets`. The items of this dictionary
are `filenames`, `raoffset` and `decoffset`. The IFU cube building code expects this dictionary to hold the information
for storing the file names and the associated ra and dec offsets. The file names should not contain the directory path.

It is assumed there exists a list of files, ra and dec offsets that are feed to this routine. The ra and dec offsets
provided in arcseconds. The cube_building code will apply the ra offsets after multiplying by the cos(crval2), where crval2 is the
It is assumed there exists a list of files, ra and dec offsets that are feed to this method. The ra and dec offsets
provided in arcseconds. The cube building code will apply the ra offsets after dividing by cos(crval2), where crval2 is the
declination center of the IFU cube.
`num` is the number of files.
y
Below `num` is the number of files.

import asdf
offsets = {}
offsets['filename'] = []
offsets['raoffset'] = []
offsets['decoffset'] = []
for i in range(num):

offsets['filename'].append(file[i])
offsets['raoffset'].append(ra_center1[i])
offsets['decoffset'].append(dec_center1[i])
Expand Down
26 changes: 25 additions & 1 deletion jwst/cube_build/ifu_cube.py
Original file line number Diff line number Diff line change
Expand Up @@ -114,7 +114,8 @@ def __init__(self,
self.naxis3 = None
self.cdelt3_normal = None
self.rot_angle = None # rotation angle between Ra-Dec and IFU local instrument plane

self.median_dec = None

self.a_min = 0
self.a_max = 0
self.b_min = 0
Expand Down Expand Up @@ -1227,6 +1228,7 @@ def setup_ifucube_wcs(self):
# _____________________________________________________________________________
self.cdelt1 = self.spatial_size
self.cdelt2 = self.spatial_size
deg2rad = math.pi / 180.0
if self.linear_wavelength:
self.cdelt3 = self.spectral_size

Expand Down Expand Up @@ -1303,6 +1305,25 @@ def setup_ifucube_wcs(self):
log.debug(f'Working on data from {this_a}, {this_b}')
n = len(self.master_table.FileMap[self.instrument][this_a][this_b])
log.debug('number of files %d', n)

# find the median center declination if we have an offset file
if self.offsets is not None:
decs = []
for k in range(n):
input_file = self.master_table.FileMap[self.instrument][this_a][this_b][k]
input_model = datamodels.open(input_file)
spatial_box = input_model.meta.wcsinfo.s_region
s = spatial_box.split(' ')
cb1 = float(s[4])
cb2 = float(s[6])
cb3 = float(s[8])
cb4 = float(s[10])
m = (cb1 + cb2 + cb3 + cb4)/4
decs.append(m)

Check warning on line 1322 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1311-L1322

Added lines #L1311 - L1322 were not covered by tests

self.median_dec = np.nanmedian(decs)
print('Median declination ', self.median_dec)

Check warning on line 1325 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1324-L1325

Added lines #L1324 - L1325 were not covered by tests

for k in range(n):
lmin = 0.0
lmax = 0.0
Expand All @@ -1322,6 +1343,7 @@ def setup_ifucube_wcs(self):
log.info("Ra and dec offset (arc seconds) applied to file :%5.2f, %5.2f, %s",

Check warning on line 1343 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1339-L1343

Added lines #L1339 - L1343 were not covered by tests
raoffset, decoffset, filename)
raoffset = raoffset/3600.0 # convert to degress
raoffset = raoffset /np.cos(self.median_dec *deg2rad)
decoffset = decoffset/3600.0

Check warning on line 1347 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1345-L1347

Added lines #L1345 - L1347 were not covered by tests
# ________________________________________________________________________________
# Find the footprint of the image
Expand Down Expand Up @@ -1731,6 +1753,7 @@ def map_miri_pixel_to_sky(self, input_model, this_par1, subtract_background,
slice_no = None # Slice number
dwave = None
corner_coord = None
deg2rad = math.pi / 180.0

Check warning on line 1756 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1756

Added line #L1756 was not covered by tests

raoffset = 0.0
decoffset = 0.0

Check warning on line 1759 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1758-L1759

Added lines #L1758 - L1759 were not covered by tests
Expand All @@ -1743,6 +1766,7 @@ def map_miri_pixel_to_sky(self, input_model, this_par1, subtract_background,
log.info("Ra and dec offset (arc seconds) applied to file :%5.2f, %5.2f, %s",

Check warning on line 1766 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1761-L1766

Added lines #L1761 - L1766 were not covered by tests
raoffset, decoffset, filename)
raoffset = raoffset/3600.0 # convert to degress
raoffset = raoffset /np.cos(self.median_dec *deg2rad)
decoffset = decoffset/3600.0

Check warning on line 1770 in jwst/cube_build/ifu_cube.py

View check run for this annotation

Codecov / codecov/patch

jwst/cube_build/ifu_cube.py#L1768-L1770

Added lines #L1768 - L1770 were not covered by tests

# check if background sky matching as been done in mrs_imatch step
Expand Down

0 comments on commit 1c5a745

Please sign in to comment.