Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

import numpy as np 

import astropy.wcs 

import astropy.units 

 

__all__ = ["WavelengthArray"] 

 

 

class WavelengthArray(np.ndarray): 

"""An array of wavelengths 

 

This subclass of `numpy.ndarray` keeps track of the construction parameters 

and uses them to create a FITS header with WCS. Specifying the wavelength 

in this way allows persisting the wavelength as a small handful of header 

keywords instead of thousands of values in a predictable series. 

 

This is a functional numpy array, although it is read-only to prevent the 

construction parameters from getting out of sync with the array values. If 

you find the need to change the values, ``copy()`` this array first. 

 

Parameters 

---------- 

minWavelength : `float` 

Minimum wavelength (nm). 

maxWavelength : `float` 

Maximum wavelength (nm). 

length : `int` 

Number of values. 

dtype : `numpy.dtype`, optional 

Data type. 

""" 

def __new__(cls, minWavelength, maxWavelength, length, dtype=np.float32): 

obj = np.linspace(minWavelength, maxWavelength, length, dtype=dtype).view(cls) 

obj.minWavelength = minWavelength 

obj.maxWavelength = maxWavelength 

obj.flags.writeable = False 

return obj 

 

def __array_finalize__(self, obj): 

"""Numpy mechanics to create array""" 

if obj is None: 

return 

self.minWavelength = getattr(obj, "minWavelength", None) 

self.maxWavelength = getattr(obj, "maxWavelength", None) 

 

def __repr__(self): 

"""String representation""" 

return (f"{type(self).__name__}({self.minWavelength}, {self.maxWavelength}, " 

f"{len(self)}, {self.dtype})") 

 

def toFitsHeader(self): 

"""Convert to a FITS header 

 

Returns 

------- 

header : `astropy.io.fits.Header` 

FITS header with WCS specifying wavelength array. 

""" 

wcs = astropy.wcs.WCS(naxis=1) 

wcs.wcs.ctype = ["WAVE"] 

wcs.wcs.cname = ["Wavelength"] 

wcs.wcs.crval = [0.5*(self.minWavelength + self.maxWavelength)] 

wcs.wcs.crpix = [0.5*len(self) + 1.5] # FITS is unit-indexed 

wcs.wcs.cunit = [astropy.units.nm] 

dWavelength = (self.maxWavelength - self.minWavelength)/(len(self) - 1) 

wcs.wcs.cdelt = [dWavelength] 

return wcs.to_header() 

 

@classmethod 

def fromFitsHeader(cls, header, length, dtype=np.float32): 

"""Construct from a FITS header 

 

Parameters 

---------- 

header : `astropy.io.fits.Header` 

FITS header with WCS specifying wavelength array. 

length : `int` 

Length of the array. 

dtype : `numpy.dtype` 

Array data type. 

 

Returns 

------- 

self : cls 

Constructed wavelength array. 

""" 

wcs = astropy.wcs.WCS(header=header) 

if len(wcs.wcs.ctype) != 1 or wcs.wcs.ctype[0] != "WAVE": 

raise RuntimeError(f"Unexpected CTYPE in header: {wcs.wcs.ctype}") 

minWavelength = wcs.pixel_to_world(1.0).to(astropy.units.nm).value 

maxWavelength = wcs.pixel_to_world(length).to(astropy.units.nm).value 

return cls(minWavelength, maxWavelength, length, dtype=dtype)