Skip to content

Commit 03482ed

Browse files
committed
Merge pull request #11 from desihub/focal_plane_transforms
Focal plane transforms
2 parents 4ff2181 + 4b23e81 commit 03482ed

12 files changed

+532
-27
lines changed

.gitignore

+1
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,7 @@ __pycache__
1414
*/cython_version.py
1515
htmlcov
1616
.coverage
17+
.cache
1718
MANIFEST
1819

1920
# Sphinx

.travis.yml

+1-1
Original file line numberDiff line numberDiff line change
@@ -110,7 +110,7 @@ install:
110110
- if [[ $SETUP_CMD == build_sphinx* ]]; then $CONDA_INSTALL numpy=$NUMPY_VERSION Sphinx matplotlib; fi
111111

112112
# COVERAGE DEPENDENCIES
113-
- if [[ $SETUP_CMD == 'test --coverage' ]]; then $PIP_INSTALL coverage coveralls; fi
113+
- if [[ $SETUP_CMD == 'test --coverage' ]]; then $PIP_INSTALL coverage==3.7.1 coveralls; fi
114114

115115
script:
116116
- python setup.py $SETUP_CMD

docs/api.rst

+7-1
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,10 @@
11
API Reference
22
=============
33

4-
.. automodapi:: specsim
4+
.. automodapi:: specsim.atmosphere
5+
.. automodapi:: specsim.driver
6+
.. automodapi:: specsim.instrument
7+
.. automodapi:: specsim.quick
8+
.. automodapi:: specsim.sources
9+
.. automodapi:: specsim.spectrum
10+
.. automodapi:: specsim.transform

specsim/__init__.py

+6-4
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,9 @@
1212

1313
# For egg_info test builds to pass, put package imports here.
1414
if not _ASTROPY_SETUP_:
15-
from .spectrum import WavelengthFunction,SpectralFluxDensity
16-
from .instrument import Instrument
17-
from .atmosphere import Atmosphere
18-
from .quick import Quick
15+
import atmosphere
16+
import instrument
17+
import quick
18+
import sources
19+
import spectrum
20+
import transform

specsim/atmosphere.py

+6-3
Original file line numberDiff line numberDiff line change
@@ -8,7 +8,8 @@
88
import numpy as np
99
import scipy.interpolate
1010

11-
from . import WavelengthFunction,SpectralFluxDensity
11+
import specsim.spectrum
12+
1213

1314
class Atmosphere(object):
1415
def __init__(self,skySpectrumFilename=None,zenithExtinctionFilename=None,
@@ -35,6 +36,8 @@ def __init__(self,skySpectrumFilename=None,zenithExtinctionFilename=None,
3536
zenithExtinctionFilename = os.path.join(basePath,
3637
'data','spectra','ZenithExtinction-KPNO.dat')
3738
# Load the tabulated sky spectrum.
38-
self.skySpectrum = SpectralFluxDensity.loadFromTextFile(skySpectrumFilename)
39+
self.skySpectrum =\
40+
specsim.spectrum.SpectralFluxDensity.loadFromTextFile(skySpectrumFilename)
3941
# Load the tabulated zenith extinction coefficients.
40-
self.zenithExtinction = WavelengthFunction.loadFromTextFile(zenithExtinctionFilename)
42+
self.zenithExtinction =\
43+
specsim.spectrum.WavelengthFunction.loadFromTextFile(zenithExtinctionFilename)

specsim/driver.py

+5-4
Original file line numberDiff line numberDiff line change
@@ -20,7 +20,8 @@
2020

2121
from astropy.utils.compat import argparse
2222

23-
import specsim as sim
23+
import specsim
24+
2425

2526
# This is a setup.py entry-point, not a standalone script.
2627
# See http://astropy.readthedocs.org/en/latest/development/scripts.html
@@ -92,7 +93,7 @@ def main(args=None):
9293
if not os.path.isfile(args.infile):
9394
print('No such infile: %s' % args.infile)
9495
return -1
95-
srcSpectrum = sim.SpectralFluxDensity.loadFromTextFile(args.infile,
96+
srcSpectrum = specsim.spectrum.SpectralFluxDensity.loadFromTextFile(args.infile,
9697
wavelengthColumn=args.infile_wavecol,valuesColumn=args.infile_fluxcol,
9798
extrapolatedValue=(0. if args.truncated else None))
9899

@@ -127,10 +128,10 @@ def main(args=None):
127128
print(specSummary)
128129

129130
# Create the default atmosphere for the requested sky conditions.
130-
atmosphere = sim.Atmosphere(skyConditions=args.sky,basePath=os.environ['SPECSIM_MODEL'])
131+
atmosphere = specsim.atmosphere.Atmosphere(skyConditions=args.sky,basePath=os.environ['SPECSIM_MODEL'])
131132

132133
# Create a quick simulator using the default instrument model.
133-
qsim = sim.Quick(atmosphere=atmosphere,basePath=os.environ['SPECSIM_MODEL'])
134+
qsim = specsim.quick.Quick(atmosphere=atmosphere,basePath=os.environ['SPECSIM_MODEL'])
134135

135136
# Initialize the simulation wavelength grid to use.
136137
if args.verbose:

specsim/instrument.py

+16-9
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,8 @@
1212
from astropy.io import fits
1313
import scipy.interpolate
1414

15-
from . import WavelengthFunction,SpectralFluxDensity
15+
import specsim.spectrum
16+
1617

1718
class Instrument(object):
1819
"""
@@ -48,16 +49,17 @@ def __init__(self,modelFilename=None,throughputPath=None,psfFile=None,basePath='
4849
# from the FITS table stored as HDU[1]. Values outside of the tabulated
4950
# wavelength range will silently extrapolate to zero.
5051
table = hduList[1].data
51-
self.throughput[camera] = WavelengthFunction(table['wavelength'],table['throughput'],
52-
extrapolatedValue=0.)
52+
self.throughput[camera] = specsim.spectrum.WavelengthFunction(
53+
table['wavelength'],table['throughput'], extrapolatedValue=0.)
5354
hduList.close()
5455
# Loop over fiberloss models present in the throughput directory.
5556
self.fiberloss = { }
5657
for fiberlossFile in glob.iglob(os.path.join(throughputPath,'fiberloss-*.dat')):
5758
# Extract the model name from the filename.
5859
model = os.path.basename(fiberlossFile)[10:-4]
5960
# Load the data from this file.
60-
self.fiberloss[model] = WavelengthFunction.loadFromTextFile(fiberlossFile,extrapolatedValue=0.)
61+
self.fiberloss[model] = specsim.spectrum.WavelengthFunction.loadFromTextFile(
62+
fiberlossFile,extrapolatedValue=0.)
6163
# Open the PSF parameter file.
6264
hduList = fits.open(psfFile)
6365
# Loop over camera bands to build linear interpolations of the PSF FWHM in the
@@ -79,13 +81,18 @@ def __init__(self,modelFilename=None,throughputPath=None,psfFile=None,basePath='
7981
table = hduList[key].data
8082
wave = table['wavelength']
8183
# Load tabulated PSF functions of wavelength.
82-
self.angstromsPerRow[band] = WavelengthFunction(wave,table['angstroms_per_row'],extrapolatedValue=0.)
83-
self.psfFWHMWavelength[band] = WavelengthFunction(wave,table['fwhm_wave'],extrapolatedValue=0.)
84-
self.psfFWHMSpatial[band] = WavelengthFunction(wave,table['fwhm_wave'],extrapolatedValue=0.)
85-
self.psfNPixelsSpatial[band] = WavelengthFunction(wave,table['neff_spatial'],extrapolatedValue=0.)
84+
self.angstromsPerRow[band] = specsim.spectrum.WavelengthFunction(
85+
wave,table['angstroms_per_row'],extrapolatedValue=0.)
86+
self.psfFWHMWavelength[band] = specsim.spectrum.WavelengthFunction(
87+
wave,table['fwhm_wave'],extrapolatedValue=0.)
88+
self.psfFWHMSpatial[band] = specsim.spectrum.WavelengthFunction(
89+
wave,table['fwhm_wave'],extrapolatedValue=0.)
90+
self.psfNPixelsSpatial[band] = specsim.spectrum.WavelengthFunction(
91+
wave,table['neff_spatial'],extrapolatedValue=0.)
8692
# Get the wavelength limits for the camera from the FITS header.
8793
waveMin,waveMax = hduList[key].header['WMIN_ALL'],hduList[key].header['WMAX_ALL']
88-
assert waveMin == wave[0] and waveMax == wave[-1], ("Inconsistent wavelength limits for %s" % key)
94+
assert waveMin == wave[0] and waveMax == wave[-1], (
95+
"Inconsistent wavelength limits for %s" % key)
8996
self.cameraWavelengthRanges.append((waveMin,waveMax))
9097
cameraMidpt.append(0.5*(waveMin+waveMax))
9198
hduList.close()

specsim/quick.py

+10-5
Original file line numberDiff line numberDiff line change
@@ -12,7 +12,10 @@
1212
from astropy import constants as const
1313
from astropy import units
1414

15-
from . import Atmosphere,Instrument,SpectralFluxDensity
15+
import specsim.spectrum
16+
import specsim.atmosphere
17+
import specsim.instrument
18+
1619

1720
class QuickCamera(object):
1821
"""
@@ -134,8 +137,10 @@ class Quick(object):
134137
using the specified base path.
135138
"""
136139
def __init__(self,atmosphere=None,instrument=None,basePath=''):
137-
self.atmosphere = atmosphere if atmosphere else Atmosphere(basePath=basePath)
138-
self.instrument = instrument if instrument else Instrument(basePath=basePath)
140+
self.atmosphere = (atmosphere if atmosphere else
141+
specsim.atmosphere.Atmosphere(basePath=basePath))
142+
self.instrument = (instrument if instrument else
143+
specsim.instrument.Instrument(basePath=basePath))
139144

140145
# Precompute the physical constant h*c in units of erg*Ang.
141146
self.hc = const.h.to(units.erg*units.s)*const.c.to(units.angstrom/units.s)
@@ -287,8 +292,8 @@ def simulate(self,sourceType,sourceSpectrum,airmass=1.,nread=1.,expTime=None,dow
287292

288293
# Convert the source to a SpectralFluxDensity if necessary, setting the flux to zero
289294
# outside the input source spectrum's range.
290-
if not isinstance(sourceSpectrum,SpectralFluxDensity):
291-
sourceSpectrum = SpectralFluxDensity(
295+
if not isinstance(sourceSpectrum,specsim.spectrum.SpectralFluxDensity):
296+
sourceSpectrum = specsim.spectrum.SpectralFluxDensity(
292297
sourceSpectrum[0],sourceSpectrum[1],extrapolatedValue=0.)
293298

294299
# Resample the source spectrum to our simulation grid, if necessary.

specsim/sources.py

+1
Original file line numberDiff line numberDiff line change
@@ -10,6 +10,7 @@
1010

1111
root2pi = np.sqrt(2*math.pi)
1212

13+
1314
def getOIIDoubletFlux(totalFlux,wave):
1415
"""
1516
Returns a vector of [OII] doublet fluxes in erg/(s*cm^2*Ang) corresponding to the

specsim/spectrum.py

+1
Original file line numberDiff line numberDiff line change
@@ -12,6 +12,7 @@
1212
from astropy import constants as const
1313
from astropy import units
1414

15+
1516
class WavelengthFunction(object):
1617
"""
1718
Represents an arbitrary function of wavelength.

specsim/tests/test_transform.py

+116
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,116 @@
1+
# Licensed under a 3-clause BSD style license - see LICENSE.rst
2+
3+
from astropy.tests.helper import pytest
4+
from ..transform import altaz_to_focalplane, focalplane_to_altaz, \
5+
observatories, sky_to_altaz, adjust_time_to_hour_angle
6+
7+
import numpy as np
8+
from astropy.time import Time
9+
from astropy.coordinates import AltAz
10+
import astropy.units as u
11+
12+
13+
def test_origin_to_focalplane():
14+
alt, az = 0.5 * u.rad, 1.5 * u.rad
15+
x, y = altaz_to_focalplane(alt, az, alt, az)
16+
assert np.allclose(x, 0 * u.rad) and np.allclose(y, 0 * u.rad)
17+
18+
19+
def test_focalplane_units():
20+
platescale = 200 * u.mm / u.deg
21+
alt, az = 0.5 * u.rad, 1.5 * u.rad
22+
x, y = altaz_to_focalplane(alt, az, alt, az, platescale=platescale)
23+
assert x.unit == u.m and y.unit == u.m
24+
alt, az = focalplane_to_altaz(x, y, alt, az, platescale=platescale)
25+
assert alt.unit == u.rad and az.unit == u.rad
26+
27+
28+
def test_shape_to_focalplane():
29+
zero = 0. * u.rad
30+
x, y = altaz_to_focalplane(zero, zero, zero, zero)
31+
assert x.shape == y.shape and x.shape == ()
32+
33+
angle = np.linspace(-0.1, +0.1, 3) * u.rad
34+
x, y = altaz_to_focalplane(angle, zero, zero, zero)
35+
assert x.shape == y.shape and x.shape == (3,)
36+
x, y = altaz_to_focalplane(angle, angle, zero, zero)
37+
assert x.shape == y.shape and x.shape == (3,)
38+
x, y = altaz_to_focalplane(angle, zero, angle, zero)
39+
assert x.shape == y.shape and x.shape == (3,)
40+
x, y = altaz_to_focalplane(angle, zero, zero, angle)
41+
assert x.shape == y.shape and x.shape == (3,)
42+
x, y = altaz_to_focalplane(angle, angle, angle, zero)
43+
assert x.shape == y.shape and x.shape == (3,)
44+
x, y = altaz_to_focalplane(angle, angle, zero, angle)
45+
assert x.shape == y.shape and x.shape == (3,)
46+
x, y = altaz_to_focalplane(angle, zero, angle, angle)
47+
assert x.shape == y.shape and x.shape == (3,)
48+
x, y = altaz_to_focalplane(zero, angle, angle, angle)
49+
assert x.shape == y.shape and x.shape == (3,)
50+
x, y = altaz_to_focalplane(angle, angle, angle, angle)
51+
assert x.shape == y.shape and x.shape == (3,)
52+
x, y = altaz_to_focalplane(angle[:, np.newaxis], angle, zero, zero)
53+
54+
assert x.shape == y.shape and x.shape == (3, 3)
55+
try:
56+
x, y = altaz_to_focalplane(angle, angle[:, np.newaxis],
57+
angle[:, np.newaxis, np.newaxis], zero)
58+
assert x.shape == y.shape and x.shape == (3, 3, 3)
59+
x, y = altaz_to_focalplane(angle, angle[:, np.newaxis],
60+
angle[:, np.newaxis, np.newaxis],
61+
angle[:, np.newaxis, np.newaxis, np.newaxis])
62+
assert x.shape == y.shape and x.shape == (3, 3, 3, 3)
63+
x, y = altaz_to_focalplane(angle, angle[:, np.newaxis],
64+
zero, angle[:, np.newaxis, np.newaxis, np.newaxis])
65+
assert x.shape == y.shape and x.shape == (3, 1, 3, 3)
66+
except RuntimeError:
67+
# These tests fails for numpy < 1.9 because np.einsum does not
68+
# broadcast correctly in this case. For details, See
69+
# https://github.com/desihub/specsim/issues/10
70+
pass
71+
72+
73+
def test_focalplane_to_origin():
74+
alt0, az0 = 0.5 * u.rad, 1.5 * u.rad
75+
alt, az = focalplane_to_altaz(0. * u.rad, 0. * u.rad, alt0, az0)
76+
assert np.allclose(alt, alt0) and np.allclose(az, az0)
77+
78+
79+
def test_focalplane_roundtrip():
80+
alt0, az0 = 0.5 * u.rad, 1.5 * u.rad
81+
x, y = -0.01 * u.rad, +0.02 * u.rad
82+
alt, az = focalplane_to_altaz(x, y, alt0, az0)
83+
x2, y2 = altaz_to_focalplane(alt, az, alt0, az0)
84+
assert np.allclose(x, x2) and np.allclose(y, y2)
85+
alt2, az2 = focalplane_to_altaz(x2, y2, alt0, az0)
86+
assert np.allclose(alt, alt2) and np.allclose(az, az2)
87+
88+
89+
def test_altaz_null():
90+
alt, az = 0.5 * u.rad, 1.5 * u.rad
91+
where = observatories['APO']
92+
when = Time(56383, format='mjd')
93+
wlen = 5400 * u.Angstrom
94+
temperature = 5 * u.deg_C
95+
pressure = 800 * u.kPa
96+
altaz_in = AltAz(alt=alt, az=az, location=where, obstime=when,
97+
obswl=wlen, temperature=temperature, pressure=pressure)
98+
altaz_out = sky_to_altaz(altaz_in, where=where, when=when,
99+
wavelength=wlen, temperature=temperature, pressure=pressure)
100+
assert np.allclose(altaz_in.alt, altaz_out.alt)
101+
assert np.allclose(altaz_in.az, altaz_out.az)
102+
103+
104+
def test_adjust_null():
105+
ra = 45 * u.deg
106+
when = Time(56383, format='mjd', location=observatories['APO'])
107+
ha = when.sidereal_time('apparent') - ra
108+
adjusted = adjust_time_to_hour_angle(when, ra, ha)
109+
assert adjusted == when
110+
111+
112+
def test_adjust_missing_longitude():
113+
ra = 45 * u.deg
114+
when = Time(56383, format='mjd', location=None)
115+
with pytest.raises(ValueError):
116+
adjusted = adjust_time_to_hour_angle(when, ra, 0.)

0 commit comments

Comments
 (0)