Coverage for python/lsst/obs/sdss/convertasTrans.py : 78%

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# LSST Data Management System
3# Copyright 2008, 2009, 2010 LSST Corporation.
4#
5# This product includes software developed by the
6# LSST Project (http://www.lsst.org/).
7#
8# This program is free software: you can redistribute it and/or modify
9# it under the terms of the GNU General Public License as published by
10# the Free Software Foundation, either version 3 of the License, or
11# (at your option) any later version.
12#
13# This program is distributed in the hope that it will be useful,
14# but WITHOUT ANY WARRANTY; without even the implied warranty of
15# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16# GNU General Public License for more details.
17#
18# You should have received a copy of the LSST License Statement and
19# the GNU General Public License along with this program. If not,
20# see <http://www.lsstcorp.org/LegalNotices/>.
21#
22import sys
24from astropy.io import fits
25import numpy as np
27import lsst.afw.image as afwImage
28from lsst.afw.geom import makeSkyWcs
29import lsst.afw.table as afwTable
30import lsst.geom as geom
31import lsst.meas.astrom.sip as sip
33deg2rad = np.pi / 180.
34rad2deg = 180. / np.pi
37class CoordinateMapper(object):
38 # COMMENT mu nu are defined as:
39 # COMMENT r'-i' < riCut:
40 # COMMENT rowm = row+dRow0+dRow1*col+dRow2*(col^2)+dRow3*(col^3)+csRow*c
41 # COMMENT colm = col+dCol0+dCol1*col+dCol2*(col^2)+dCol3*(col^3)+csCol*c
42 # COMMENT r'-i' >= riCut
43 # COMMENT rowm = row+dRow0+dRow1*col+dRow2*(col^2)+dRow3*(col^3)+ccRow
44 # COMMENT colm = col+dCol0+dCol1*col+dCol2*(col^2)+dCol3*(col^3)+ccCol
45 # COMMENT mu = a + b * rowm + c * colm
46 # COMMENT nu = d + e * rowm + f * colm
48 # We will do an "average" mapping for an r-i=0 object so that the
49 # last cs/cc terms are not necessary
51 def __init__(self, node_rad, incl_rad, dRow0, dRow1, dRow2, dRow3, dCol0, dCol1, dCol2, dCol3,
52 a, b, c, d, e, f, cOffset=+0.5):
53 # Here cOffset reflects the differences between SDSS coords
54 # (LLC = 0.5,0.5) and LSST coords (LLC = 0,0). If SDSS
55 # measures an object centered at (0.5,0.5) then LSST will
56 # measure it at coordinate (0,0), but when using SDSS
57 # astrometry we need to evaluate the equations at (x+0.5,
58 # y+0.5)
60 self.node_rad = node_rad
61 self.incl_rad = incl_rad
63 self.dRow0 = dRow0
64 self.dRow1 = dRow1
65 self.dRow2 = dRow2
66 self.dRow3 = dRow3
68 self.dCol0 = dCol0
69 self.dCol1 = dCol1
70 self.dCol2 = dCol2
71 self.dCol3 = dCol3
73 self.a = a
74 self.b = b
75 self.c = c
76 self.d = d
77 self.e = e
78 self.f = f
80 self.cOff = cOffset
82 def xyToMuNu(self, x, y):
83 rowm = (y+self.cOff) + self.dRow0 + self.dRow1*(x+self.cOff) + self.dRow2*((x+self.cOff)**2) + \
84 self.dRow3*((x+self.cOff)**3)
85 colm = (x+self.cOff) + self.dCol0 + self.dCol1*(x+self.cOff) + self.dCol2*((x+self.cOff)**2) + \
86 self.dCol3*((x+self.cOff)**3)
88 mu_deg = self.a + self.b * rowm + self.c * colm
89 nu_deg = self.d + self.e * rowm + self.f * colm
90 mu_rad = mu_deg * deg2rad
91 nu_rad = nu_deg * deg2rad
93 return mu_rad, nu_rad
95 def muNuToRaDec(self, mu_rad, nu_rad):
96 x2 = np.cos(mu_rad - self.node_rad) * np.cos(nu_rad)
97 y2 = np.sin(mu_rad - self.node_rad) * np.cos(nu_rad)
98 z2 = np.sin(nu_rad)
99 y1 = y2 * np.cos(self.incl_rad) - z2 * np.sin(self.incl_rad)
100 z1 = y2 * np.sin(self.incl_rad) + z2 * np.cos(self.incl_rad)
102 ra_rad = np.arctan2(y1, x2) + self.node_rad
103 dec_rad = np.arcsin(z1)
105 return ra_rad, dec_rad
107 def xyToRaDec(self, x, y):
108 mu_rad, nu_rad = self.xyToMuNu(x, y)
109 return self.muNuToRaDec(mu_rad, nu_rad)
112def createWcs(x, y, mapper, order=4, cOffset=1.0):
113 # Here cOffset reflects the differences between FITS coords (LLC =
114 # 1,1) and LSST coords (LLC = 0,0). That is, when creating a Wcs
115 # from scratch, we need to evaluate our WCS at coordinate 0,0 to
116 # create CRVAL, but set CRPIX to 1,1.
118 ra_rad, dec_rad = mapper.xyToRaDec(x, y)
120 # Minimial table for sky coordinates
121 catTable = afwTable.SimpleTable.make(afwTable.SimpleTable.makeMinimalSchema())
123 # Minimial table + centroids for focal plane coordintes
124 srcSchema = afwTable.SourceTable.makeMinimalSchema()
125 centroidKey = afwTable.Point2DKey.addFields(srcSchema, "centroid", "centroid", "pixel")
127 srcTable = afwTable.SourceTable.make(srcSchema)
128 srcTable.defineCentroid("centroid")
130 matches = []
131 for i in range(len(x)):
132 src = srcTable.makeRecord()
133 src.set(centroidKey.getX(), x[i])
134 src.set(centroidKey.getY(), y[i])
136 cat = catTable.makeRecord()
137 cat.set(catTable.getCoordKey().getRa(), geom.Angle(ra_rad[i], geom.radians))
138 cat.set(catTable.getCoordKey().getDec(), geom.Angle(dec_rad[i], geom.radians))
140 mat = afwTable.ReferenceMatch(cat, src, 0.0)
141 matches.append(mat)
143 # Need to make linear Wcs around which to expand solution
145 # CRPIX1 = Column Pixel Coordinate of Ref. Pixel
146 # CRPIX2 = Row Pixel Coordinate of Ref. Pixel
147 crpix = geom.Point2D(x[0] + cOffset, y[0] + cOffset)
149 # CRVAL1 = RA at Reference Pixel
150 # CRVAL2 = DEC at Reference Pixel
151 crval = geom.SpherePoint(ra_rad[0], dec_rad[0], geom.radians)
153 # CD1_1 = RA degrees per column pixel
154 # CD1_2 = RA degrees per row pixel
155 # CD2_1 = DEC degrees per column pixel
156 # CD2_2 = DEC degrees per row pixel
157 LLl = mapper.xyToRaDec(0., 0.)
158 ULl = mapper.xyToRaDec(0., 1.)
159 LRl = mapper.xyToRaDec(1., 0.)
161 LLc = geom.SpherePoint(LLl[0], LLl[1], geom.radians)
162 ULc = geom.SpherePoint(ULl[0], ULl[1], geom.radians)
163 LRc = geom.SpherePoint(LRl[0], LRl[1], geom.radians)
165 cdN_1 = LLc.getTangentPlaneOffset(LRc)
166 cdN_2 = LLc.getTangentPlaneOffset(ULc)
167 cd1_1, cd2_1 = cdN_1[0].asDegrees(), cdN_1[1].asDegrees()
168 cd1_2, cd2_2 = cdN_2[0].asDegrees(), cdN_2[1].asDegrees()
170 cdMatrix = np.array([cd1_1, cd2_1, cd1_2, cd2_2], dtype=float)
171 cdMatrix.shape = (2, 2)
173 linearWcs = makeSkyWcs(crpix=crpix, crval=crval, cdMatrix=cdMatrix)
174 wcs = sip.makeCreateWcsWithSip(matches, linearWcs, order).getNewWcs()
176 return wcs
179def validate(xs, ys, mapper, wcs):
180 dists = []
181 for i in range(len(xs)):
182 tuple1 = mapper.xyToRaDec(xs[i], ys[i])
183 coord1 = geom.SpherePoint(tuple1[0], tuple1[1], geom.radians)
184 coord2 = wcs.pixelToSky(xs[i], ys[i])
185 dist = coord1.separation(coord2).asArcseconds()
186 dists.append(dist)
188 print(np.mean(dists), np.std(dists))
191def convertasTrans(infile, filt, camcol, field, stepSize=50, doValidate=False):
192 with fits.open(infile) as hdulist:
193 t0 = hdulist[0].header['ccdarray']
194 if t0 != 'photo': 194 ↛ 195line 194 didn't jump to line 195, because the condition on line 194 was never true
195 raise RuntimeError('*** Cannot support ccdarray: %s' % (t0,))
197 camcols = hdulist[0].header['camcols']
198 filters = hdulist[0].header['filters']
199 node_deg = hdulist[0].header['node']
200 incl_deg = hdulist[0].header['incl']
201 node_rad = node_deg * deg2rad
202 incl_rad = incl_deg * deg2rad
204 cList = [int(cc) for cc in camcols.split()]
205 fList = filters.split()
207 try:
208 cIdx = cList.index(camcol)
209 except Exception:
210 print("Cannot extract data for camcol %s" % (camcol))
211 return None
213 try:
214 fIdx = fList.index(filt)
215 except Exception:
216 print("Cannot extract data for filter %s" % (filt))
217 return None
219 ext = cIdx * len(fList) + (fIdx + 1)
220 ehdr = hdulist[ext].header
221 edat = hdulist[ext].data
223 if ehdr['CAMCOL'] != camcol or ehdr['FILTER'] != filt: 223 ↛ 224line 223 didn't jump to line 224, because the condition on line 223 was never true
224 print("Extracted incorrect header; fix me")
225 return None
227 fields = edat.field('field').tolist()
228 try:
229 fIdx = fields.index(field)
230 except Exception:
231 print("Cannot extract data for field %d" % (field))
232 return None
234 dRow0 = edat.field('dRow0')[fIdx]
235 dRow1 = edat.field('dRow1')[fIdx]
236 dRow2 = edat.field('dRow2')[fIdx]
237 dRow3 = edat.field('dRow3')[fIdx]
239 dCol0 = edat.field('dCol0')[fIdx]
240 dCol1 = edat.field('dCol1')[fIdx]
241 dCol2 = edat.field('dCol2')[fIdx]
242 dCol3 = edat.field('dCol3')[fIdx]
244 a = edat.field('a')[fIdx]
245 b = edat.field('b')[fIdx]
246 c = edat.field('c')[fIdx]
247 d = edat.field('d')[fIdx]
248 e = edat.field('e')[fIdx]
249 f = edat.field('f')[fIdx]
251 # We need to fit for a TAN-SIP
252 x = np.arange(0, 1489+stepSize, stepSize)
253 y = np.arange(0, 2048+stepSize, stepSize)
254 coords = np.meshgrid(x, y)
255 xs = np.ravel(coords[0]).astype(np.float)
256 ys = np.ravel(coords[1]).astype(np.float)
257 mapper = CoordinateMapper(node_rad, incl_rad, dRow0, dRow1, dRow2, dRow3, dCol0, dCol1, dCol2, dCol3,
258 a, b, c, d, e, f)
259 wcs = createWcs(xs, ys, mapper)
261 if doValidate: 261 ↛ 262line 261 didn't jump to line 262, because the condition on line 261 was never true
262 validate(xs, ys, mapper, wcs)
264 return wcs
267if __name__ == '__main__': 267 ↛ 268line 267 didn't jump to line 268, because the condition on line 267 was never true
268 infile = sys.argv[1]
269 filt = sys.argv[2]
270 camcol = int(sys.argv[3])
271 field = int(sys.argv[4])
273 wcs = convertasTrans(infile, filt, camcol, field, doValidate=True)
275 if len(sys.argv) > 5:
276 fpC = sys.argv[5]
277 image = afwImage.ImageF(fpC)
278 mi = afwImage.MaskedImageF(image)
279 exp = afwImage.ExposureF(mi, wcs)
280 exp.writeFits("/tmp/exp.fits")