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""" 

2This file contains coordinate transformation methods that are very thin wrappers 

3of palpy methods, or that have no dependence on palpy at all 

4""" 

5from __future__ import division 

6 

7import numpy as np 

8import numbers 

9import palpy 

10 

11from lsst.sims.utils.CodeUtilities import _validate_inputs 

12 

13__all__ = ["_galacticFromEquatorial", "galacticFromEquatorial", 

14 "_equatorialFromGalactic", "equatorialFromGalactic", 

15 "sphericalFromCartesian", "cartesianFromSpherical", 

16 "xyz_from_ra_dec", "_xyz_from_ra_dec", 

17 "_ra_dec_from_xyz", "ra_dec_from_xyz", 

18 "xyz_angular_radius", "_xyz_angular_radius", 

19 "rotationMatrixFromVectors", 

20 "rotAboutZ", "rotAboutY", "rotAboutX", 

21 "equationOfEquinoxes", "calcGmstGast", "calcLmstLast", 

22 "angularSeparation", "_angularSeparation", "haversine", 

23 "arcsecFromRadians", "radiansFromArcsec", 

24 "arcsecFromDegrees", "degreesFromArcsec"] 

25 

26 

27def calcLmstLast(mjd, longRad): 

28 """ 

29 calculates local mean sidereal time and local apparent sidereal time 

30 

31 @param [in] mjd is the universal time (UT1) expressed as an MJD. 

32 This can be a numpy array or a single value. 

33 

34 @param [in] longRad is the longitude in radians (positive east of the prime meridian) 

35 This can be numpy array or a single value. If a numpy array, should have the same length as mjd. In that 

36 case, each longRad will be applied only to the corresponding mjd. 

37 

38 @param [out] lmst is the local mean sidereal time in hours 

39 

40 @param [out] last is the local apparent sideral time in hours 

41 """ 

42 mjdIsArray = False 

43 longRadIsArray = False 

44 if isinstance(mjd, np.ndarray): 

45 mjdIsArray = True 

46 

47 if isinstance(longRad, np.ndarray): 

48 longRadIsArray = True 

49 

50 if longRadIsArray and mjdIsArray: 

51 if len(longRad) != len(mjd): 

52 raise RuntimeError("In calcLmstLast mjd and longRad have different lengths") 

53 

54 valid_type = False 

55 if isinstance(mjd, np.ndarray) and isinstance(longRad, np.ndarray): 

56 valid_type = True 

57 elif isinstance(mjd, np.ndarray) and isinstance(longRad, numbers.Number): 

58 valid_type = True 

59 elif isinstance(mjd, numbers.Number) and isinstance(longRad, numbers.Number): 

60 valid_type = True 

61 

62 if not valid_type: 

63 msg = "Valid input types for calcLmstLast are:\n" \ 

64 "mjd and longRad as numpy arrays of the same length\n" \ 

65 "mjd as a numpy array and longRad as a number\n" \ 

66 "mjd as a number and longRad as a number\n" \ 

67 "You gave mjd: %s\n" % type(mjd) \ 

68 + "and longRad: %s\n" % type(longRad) 

69 

70 raise RuntimeError(msg) 

71 

72 longDeg0 = np.degrees(longRad) 

73 longDeg0 %= 360.0 

74 

75 if longRadIsArray: 

76 longDeg = np.where(longDeg0 > 180.0, longDeg0 - 360.0, longDeg0) 

77 else: 

78 if longDeg0 > 180.: 

79 longDeg = longDeg0 - 360. 

80 else: 

81 longDeg = longDeg0 

82 

83 hrs = longDeg / 15.0 

84 gmstgast = calcGmstGast(mjd) 

85 lmst = gmstgast[0] + hrs 

86 last = gmstgast[1] + hrs 

87 lmst %= 24. 

88 last %= 24. 

89 return lmst, last 

90 

91 

92def galacticFromEquatorial(ra, dec): 

93 '''Convert RA,Dec (J2000) to Galactic Coordinates 

94 

95 @param [in] ra is right ascension in degrees, either a number or a numpy array 

96 

97 @param [in] dec is declination in degrees, either a number or a numpy array 

98 

99 @param [out] gLong is galactic longitude in degrees 

100 

101 @param [out] gLat is galactic latitude in degrees 

102 ''' 

103 

104 gLong, gLat = _galacticFromEquatorial(np.radians(ra), np.radians(dec)) 

105 return np.degrees(gLong), np.degrees(gLat) 

106 

107 

108def _galacticFromEquatorial(ra, dec): 

109 '''Convert RA,Dec (J2000) to Galactic Coordinates 

110 

111 All angles are in radians 

112 

113 @param [in] ra is right ascension in radians, either a number or a numpy array 

114 

115 @param [in] dec is declination in radians, either a number or a numpy array 

116 

117 @param [out] gLong is galactic longitude in radians 

118 

119 @param [out] gLat is galactic latitude in radians 

120 ''' 

121 

122 if isinstance(ra, np.ndarray): 

123 gLong, gLat = palpy.eqgalVector(ra, dec) 

124 else: 

125 gLong, gLat = palpy.eqgal(ra, dec) 

126 

127 return gLong, gLat 

128 

129 

130def equatorialFromGalactic(gLong, gLat): 

131 '''Convert Galactic Coordinates to RA, dec (J2000) 

132 

133 @param [in] gLong is galactic longitude in degrees, either a number or a numpy array 

134 (0 <= gLong <= 360.) 

135 

136 @param [in] gLat is galactic latitude in degrees, either a number or a numpy array 

137 (-90. <= gLat <= 90.) 

138 

139 @param [out] ra is right ascension in degrees 

140 

141 @param [out] dec is declination in degrees 

142 ''' 

143 

144 ra, dec = _equatorialFromGalactic(np.radians(gLong), np.radians(gLat)) 

145 return np.degrees(ra), np.degrees(dec) 

146 

147 

148def _equatorialFromGalactic(gLong, gLat): 

149 '''Convert Galactic Coordinates to RA, dec (J2000) 

150 

151 @param [in] gLong is galactic longitude in radians, either a number or a numpy array 

152 (0 <= gLong <= 2*pi) 

153 

154 @param [in] gLat is galactic latitude in radians, either a number or a numpy array 

155 (-pi/2 <= gLat <= pi/2) 

156 

157 @param [out] ra is right ascension in radians (J2000) 

158 

159 @param [out] dec is declination in radians (J2000) 

160 ''' 

161 

162 if isinstance(gLong, np.ndarray): 

163 ra, dec = palpy.galeqVector(gLong, gLat) 

164 else: 

165 ra, dec = palpy.galeq(gLong, gLat) 

166 

167 return ra, dec 

168 

169 

170def cartesianFromSpherical(longitude, latitude): 

171 """ 

172 Transforms between spherical and Cartesian coordinates. 

173 

174 @param [in] longitude is a numpy array or a number in radians 

175 

176 @param [in] latitude is a numpy array or number in radians 

177 

178 @param [out] a numpy array of the (three-dimensional) cartesian coordinates on a unit sphere. 

179 

180 if inputs are numpy arrays: 

181 output[i][0] will be the x-coordinate of the ith point 

182 output[i][1] will be the y-coordinate of the ith point 

183 output[i][2] will be the z-coordinate of the ith point 

184 

185 All angles are in radians 

186 

187 Also, look at xyz_from_ra_dec(). 

188 """ 

189 return _xyz_from_ra_dec(longitude, latitude).transpose() 

190 

191 

192def sphericalFromCartesian(xyz): 

193 """ 

194 Transforms between Cartesian and spherical coordinates 

195 

196 @param [in] xyz is a numpy array of points in 3-D space. 

197 Each row is a different point. 

198 

199 @param [out] returns longitude and latitude 

200 

201 All angles are in radians 

202 

203 Also, look at ra_dec_from_xyz(). 

204 """ 

205 if not isinstance(xyz, np.ndarray): 

206 raise RuntimeError("You need to pass a numpy array to sphericalFromCartesian") 

207 

208 if len(xyz.shape) > 1: 

209 return _ra_dec_from_xyz(xyz[:, 0], xyz[:, 1], xyz[:, 2]) 

210 else: 

211 return _ra_dec_from_xyz(xyz[0], xyz[1], xyz[2]) 

212 

213 

214def xyz_from_ra_dec(ra, dec): 

215 """ 

216 Utility to convert RA,dec positions in x,y,z space. 

217 

218 Parameters 

219 ---------- 

220 ra : float or array 

221 RA in degrees 

222 dec : float or array 

223 Dec in degrees 

224 

225 Returns 

226 ------- 

227 x,y,z : floats or arrays 

228 The position of the given points on the unit sphere. 

229 """ 

230 return _xyz_from_ra_dec(np.radians(ra), np.radians(dec)) 

231 

232 

233def _xyz_from_ra_dec(ra, dec): 

234 """ 

235 Utility to convert RA,dec positions in x,y,z space. 

236 

237 Parameters 

238 ---------- 

239 ra : float or array 

240 RA in radians 

241 dec : float or array 

242 Dec in radians 

243 

244 Returns 

245 ------- 

246 x,y,z : floats or arrays 

247 The position of the given points on the unit sphere. 

248 """ 

249 # It is ok to mix floats and numpy arrays. 

250 

251 cosDec = np.cos(dec) 

252 return np.array([np.cos(ra) * cosDec, np.sin(ra) * cosDec, np.sin(dec)]) 

253 

254 

255def _ra_dec_from_xyz(x, y, z): 

256 """ 

257 Utility to convert x,y,z Cartesian coordinates to RA, dec positions in space. 

258 

259 Parameters 

260 ---------- 

261 x : float or array 

262 The position on the x-axis of the given points on the unit sphere 

263 y : float or array 

264 The position on the y-axis of the given points on the unit sphere 

265 z : float or array 

266 The position on the z-axis of the given points on the unit sphere 

267 

268 Returns 

269 ------- 

270 ra, dec : floats or arrays 

271 Ra and dec coordinates in radians. 

272 """ 

273 rad = np.sqrt(x**2 + y**2 + z**2) 

274 ra = np.arctan2(y, x) 

275 dec = np.arcsin(z / rad) 

276 

277 return ra, dec 

278 

279 

280def ra_dec_from_xyz(x, y, z): 

281 """ 

282 Utility to convert x,y,z Cartesian coordinates to RA, dec positions in space. 

283 

284 Parameters 

285 ---------- 

286 x : float or array 

287 The position on the x-axis of the given points on the unit sphere 

288 y : float or array 

289 The position on the y-axis of the given points on the unit sphere 

290 z : float or array 

291 The position on the z-axis of the given points on the unit sphere 

292 

293 Returns 

294 ------- 

295 ra, dec : floats or arrays 

296 Ra and dec coordinates in degrees. 

297 """ 

298 

299 return np.degrees(_ra_dec_from_xyz(x, y, z)) 

300 

301 

302def xyz_angular_radius(radius=1.75): 

303 """ 

304 Convert an angular radius into a physical radius for a kdtree search. 

305 

306 Parameters 

307 ---------- 

308 radius : float 

309 Radius in degrees. 

310 

311 Returns 

312 ------- 

313 radius : float 

314 """ 

315 return _xyz_angular_radius(np.radians(radius)) 

316 

317 

318def _xyz_angular_radius(radius): 

319 """ 

320 Convert an angular radius into a physical radius for a kdtree search. 

321 

322 Parameters 

323 ---------- 

324 radius : float 

325 Radius in radians. 

326 

327 Returns 

328 ------- 

329 radius : float 

330 """ 

331 x0, y0, z0 = (1, 0, 0) 

332 x1, y1, z1 = _xyz_from_ra_dec(radius, 0) 

333 result = np.sqrt((x1-x0)**2+(y1-y0)**2+(z1-z0)**2) 

334 return result 

335 

336 

337def zRotationMatrix(theta): 

338 cc = np.cos(theta) 

339 ss = np.sin(theta) 

340 return np.array([[cc, -ss, 0.0], 

341 [ss, cc, 0.0], 

342 [0.0, 0.0, 1.0]]) 

343 

344 

345def rotAboutZ(vec, theta): 

346 """ 

347 Rotate a Cartesian vector an angle theta about the z axis. 

348 Theta is in radians. 

349 Positive theta rotates +x towards +y. 

350 """ 

351 return np.dot(zRotationMatrix(theta), vec.transpose()).transpose() 

352 

353 

354def yRotationMatrix(theta): 

355 cc = np.cos(theta) 

356 ss = np.sin(theta) 

357 return np.array([[cc, 0.0, ss], 

358 [0.0, 1.0, 0.0], 

359 [-ss, 0.0, cc]]) 

360 

361 

362def rotAboutY(vec, theta): 

363 """ 

364 Rotate a Cartesian vector an angle theta about the y axis. 

365 Theta is in radians. 

366 Positive theta rotates +x towards -z. 

367 """ 

368 return np.dot(yRotationMatrix(theta), vec.transpose()).transpose() 

369 

370 

371def xRotationMatrix(theta): 

372 cc = np.cos(theta) 

373 ss = np.sin(theta) 

374 

375 return np.array([[1.0, 0.0, 0.0], 

376 [0.0, cc, -ss], 

377 [0.0, ss, cc]]) 

378 

379 

380def rotAboutX(vec, theta): 

381 """ 

382 Rotate a Cartesian vector an angle theta about the x axis. 

383 Theta is in radians. 

384 Positive theta rotates +y towards +z. 

385 """ 

386 return np.dot(xRotationMatrix(theta), vec.transpose()).transpose() 

387 

388 

389def rotationMatrixFromVectors(v1, v2): 

390 ''' 

391 Given two vectors v1,v2 calculate the rotation matrix for v1->v2 using the axis-angle approach 

392 

393 @param [in] v1, v2 are two Cartesian unit vectors (in three dimensions) 

394 

395 @param [out] rot is the rotation matrix that rotates from one to the other 

396 

397 ''' 

398 

399 if np.abs(np.sqrt(np.dot(v1, v1)) - 1.0) > 0.01: 

400 raise RuntimeError("v1 in rotationMatrixFromVectors is not a unit vector") 

401 

402 if np.abs(np.sqrt(np.dot(v2, v2)) - 1.0) > 0.01: 

403 raise RuntimeError("v2 in rotationMatrixFromVectors is not a unit vector") 

404 

405 # Calculate the axis of rotation by the cross product of v1 and v2 

406 cross = np.cross(v1, v2) 

407 cross = cross / np.sqrt(np.dot(cross, cross)) 

408 

409 # calculate the angle of rotation via dot product 

410 angle = np.arccos(np.dot(v1, v2)) 

411 sinDot = np.sin(angle) 

412 cosDot = np.cos(angle) 

413 

414 # calculate the corresponding rotation matrix 

415 # http://en.wikipedia.org/wiki/Rotation_matrix#Rotation_matrix_from_axis_and_angle 

416 rot = [[cosDot + cross[0] * cross[0] * (1 - cosDot), -cross[2] * sinDot + 

417 (1 - cosDot) * cross[0] * cross[1], 

418 cross[1] * sinDot + (1 - cosDot) * cross[0] * cross[2]], 

419 [cross[2] * sinDot + (1 - cosDot) * cross[0] * cross[1], cosDot + 

420 (1 - cosDot) * cross[1] * cross[1], 

421 -cross[0] * sinDot + (1 - cosDot) * cross[1] * cross[2]], 

422 [-cross[1] * sinDot + (1 - cosDot) * cross[0] * cross[2], cross[0] * sinDot + 

423 (1 - cosDot) * cross[1] * cross[2], cosDot + (1 - cosDot) * (cross[2] * cross[2])]] 

424 

425 return rot 

426 

427 

428def equationOfEquinoxes(d): 

429 """ 

430 The equation of equinoxes. See http://aa.usno.navy.mil/faq/docs/GAST.php 

431 

432 @param [in] d is either a numpy array or a number that is Terrestrial Time 

433 expressed as an MJD 

434 

435 @param [out] the equation of equinoxes in radians. 

436 """ 

437 

438 if isinstance(d, np.ndarray): 

439 return palpy.eqeqxVector(d) 

440 else: 

441 return palpy.eqeqx(d) 

442 

443 

444def calcGmstGast(mjd): 

445 """ 

446 Compute Greenwich mean sidereal time and Greenwich apparent sidereal time 

447 see: From http://aa.usno.navy.mil/faq/docs/GAST.php 

448 

449 @param [in] mjd is the universal time (UT1) expressed as an MJD 

450 

451 @param [out] gmst Greenwich mean sidereal time in hours 

452 

453 @param [out] gast Greenwich apparent sidereal time in hours 

454 """ 

455 

456 date = np.floor(mjd) 

457 ut1 = mjd - date 

458 if isinstance(mjd, np.ndarray): 

459 gmst = palpy.gmstaVector(date, ut1) 

460 else: 

461 gmst = palpy.gmsta(date, ut1) 

462 

463 eqeq = equationOfEquinoxes(mjd) 

464 gast = gmst + eqeq 

465 

466 gmst = gmst * 24.0 / (2.0 * np.pi) 

467 gmst %= 24.0 

468 

469 gast = gast * 24.0 / (2.0 * np.pi) 

470 gast %= 24.0 

471 

472 return gmst, gast 

473 

474 

475def _angularSeparation(long1, lat1, long2, lat2): 

476 """ 

477 Angular separation between two points in radians 

478 

479 Parameters 

480 ---------- 

481 long1 is the first longitudinal coordinate in radians 

482 

483 lat1 is the first latitudinal coordinate in radians 

484 

485 long2 is the second longitudinal coordinate in radians 

486 

487 lat2 is the second latitudinal coordinate in radians 

488 

489 Returns 

490 ------- 

491 The angular separation between the two points in radians 

492 

493 Calculated based on the haversine formula 

494 From http://en.wikipedia.org/wiki/Haversine_formula 

495 """ 

496 are_arrays_1 = _validate_inputs([long1, lat1], 

497 ['long1', 'lat1'], 

498 'angularSeparation') 

499 

500 are_arrays_2 = _validate_inputs([long2, lat2], 

501 ['long2', 'lat2'], 

502 'angularSeparation') 

503 

504 # The code below is necessary because the call to np.radians() in 

505 # angularSeparation() will automatically convert floats 

506 # into length 1 numpy arrays, confusing validate_inputs() 

507 if are_arrays_1 and len(long1) == 1: 

508 are_arrays_1 = False 

509 long1 = long1[0] 

510 lat1 = lat1[0] 

511 

512 if are_arrays_2 and len(long2) == 1: 

513 are_arrays_2 = False 

514 long2 = long2[0] 

515 lat2 = lat2[0] 

516 

517 if are_arrays_1 and are_arrays_2: 

518 if len(long1) != len(long2): 

519 raise RuntimeError("You need to pass arrays of the same length " 

520 "as arguments to angularSeparation()") 

521 

522 t1 = np.sin(lat2/2.0 - lat1/2.0)**2 

523 t2 = np.cos(lat1)*np.cos(lat2)*np.sin(long2/2.0 - long1/2.0)**2 

524 _sum = t1 + t2 

525 

526 if isinstance(_sum, numbers.Number): 

527 if _sum<0.0: 

528 _sum = 0.0 

529 else: 

530 _sum = np.where(_sum<0.0, 0.0, _sum) 

531 

532 return 2.0*np.arcsin(np.sqrt(_sum)) 

533 

534def angularSeparation(long1, lat1, long2, lat2): 

535 """ 

536 Angular separation between two points in degrees 

537 

538 Parameters 

539 ---------- 

540 long1 is the first longitudinal coordinate in degrees 

541 

542 lat1 is the first latitudinal coordinate in degrees 

543 

544 long2 is the second longitudinal coordinate in degrees 

545 

546 lat2 is the second latitudinal coordinate in degrees 

547 

548 Returns 

549 ------- 

550 The angular separation between the two points in degrees 

551 """ 

552 return np.degrees(_angularSeparation(np.radians(long1), 

553 np.radians(lat1), 

554 np.radians(long2), 

555 np.radians(lat2))) 

556 

557 

558def haversine(long1, lat1, long2, lat2): 

559 """ 

560 DEPRECATED; use angularSeparation() instead 

561 

562 Return the angular distance between two points in radians 

563 

564 @param [in] long1 is the longitude of point 1 in radians 

565 

566 @param [in] lat1 is the latitude of point 1 in radians 

567 

568 @param [in] long2 is the longitude of point 2 in radians 

569 

570 @param [in] lat2 is the latitude of point 2 in radians 

571 

572 @param [out] the angular separation between points 1 and 2 in radians 

573 """ 

574 return _angularSeparation(long1, lat1, long2, lat2) 

575 

576 

577def arcsecFromRadians(value): 

578 """ 

579 Convert an angle in radians to arcseconds 

580 

581 Note: if you input None, you will get None back 

582 """ 

583 if value is None: 

584 return None 

585 

586 return 3600.0 * np.degrees(value) 

587 

588 

589def radiansFromArcsec(value): 

590 """ 

591 Convert an angle in arcseconds to radians 

592 

593 Note: if you input None, you will get None back 

594 """ 

595 if value is None: 

596 return None 

597 

598 return np.radians(value / 3600.0) 

599 

600 

601def arcsecFromDegrees(value): 

602 """ 

603 Convert an angle in degrees to arcseconds 

604 

605 Note: if you input None, you will get None back 

606 """ 

607 if value is None: 

608 return None 

609 

610 return 3600.0 * value 

611 

612 

613def degreesFromArcsec(value): 

614 """ 

615 Convert an angle in arcseconds to degrees 

616 

617 Note: if you input None, you will get None back 

618 """ 

619 if value is None: 

620 return None 

621 

622 return value / 3600.0