LSSTApplications  1.1.2+25,10.0+13,10.0+132,10.0+133,10.0+224,10.0+41,10.0+8,10.0-1-g0f53050+14,10.0-1-g4b7b172+19,10.0-1-g61a5bae+98,10.0-1-g7408a83+3,10.0-1-gc1e0f5a+19,10.0-1-gdb4482e+14,10.0-11-g3947115+2,10.0-12-g8719d8b+2,10.0-15-ga3f480f+1,10.0-2-g4f67435,10.0-2-gcb4bc6c+26,10.0-28-gf7f57a9+1,10.0-3-g1bbe32c+14,10.0-3-g5b46d21,10.0-4-g027f45f+5,10.0-4-g86f66b5+2,10.0-4-gc4fccf3+24,10.0-40-g4349866+2,10.0-5-g766159b,10.0-5-gca2295e+25,10.0-6-g462a451+1
LSSTDataManagementBasePackage
Wcs.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
25 #include <iostream>
26 #include <sstream>
27 #include <cmath>
28 #include <cstring>
29 
30 #include "boost/format.hpp"
31 
32 #include "wcslib/wcs.h"
33 #include "wcslib/wcsfix.h"
34 #include "wcslib/wcshdr.h"
35 
36 #include "lsst/daf/base.h"
37 #include "lsst/daf/base/Citizen.h"
40 #include "lsst/pex/exceptions.h"
42 #include "lsst/afw/image/Wcs.h"
43 #include "lsst/afw/coord/Coord.h"
44 #include "lsst/afw/geom/Angle.h"
48 
49 namespace except = lsst::pex::exceptions;
50 namespace afwImg = lsst::afw::image;
51 namespace afwCoord = lsst::afw::coord;
52 namespace afwGeom = lsst::afw::geom;
53 
54 
55 using namespace std;
56 
63 
64 //The amount of space allocated to strings in wcslib
65 const int STRLEN = 72;
66 
67 //Set internal params for wcslib
69 {
70  _wcsfixCtrl = // ctrl for wcsfix
71  2; // Translate "H" to "h"
72  _wcshdrCtrl = // ctrl for wcspih
73  2; // Report each rejected keyrecord and the reason why it was rejected
74  _relax = // relax parameter for wcspih;
75  WCSHDR_all; // Accept all extensions recognized by the parser
76 }
77 
78 const int lsstToFitsPixels = +1;
79 const int fitsToLsstPixels = -1;
80 
81 //
82 // Constructors
83 //
84 
85 
88  daf::base::Citizen(typeid(this)),
89  _wcsInfo(NULL), _nWcsInfo(0), _relax(0), _wcsfixCtrl(0), _wcshdrCtrl(0), _nReject(0),
90  _coordSystem(static_cast<afwCoord::CoordSystem>(-1)) {
92  _initWcs();
93 }
94 
95 
99  daf::base::Citizen(typeid(this)),
100  _wcsInfo(NULL),
101  _nWcsInfo(0),
102  _relax(0),
103  _wcsfixCtrl(0),
104  _wcshdrCtrl(0),
105  _nReject(0),
106  _coordSystem(static_cast<afwCoord::CoordSystem>(-1))
107 {
109 
110  initWcsLibFromFits(fitsMetadata);
111  _initWcs();
112 }
113 
114 /*
115  * Set some internal variables that we need to refer to
116  */
118 {
119  if (_wcsInfo) {
121 
122  // tell WCSlib that values have been updated
123  _wcsInfo->flag = 0;
124  // and then tell it to do its internal magic.
125  int status = wcsset(_wcsInfo);
126  if (status != 0) {
127  throw LSST_EXCEPT(except::RuntimeError,
128  (boost::format("Failed to setup wcs structure with wcsset. Status %d: %s") %
129  status % wcs_errmsg[status] ).str());
130  }
131  }
132 }
133 
148 Wcs::Wcs(GeomPoint const & crval, GeomPoint const & crpix, Eigen::Matrix2d const & CD,
149  std::string const & ctype1, std::string const & ctype2,
150  double equinox, std::string const & raDecSys,
151  std::string const & cunits1, std::string const & cunits2
152 ):
153  daf::base::Citizen(typeid(this)),
154  _wcsInfo(NULL),
155  _nWcsInfo(0),
156  _relax(0),
157  _wcsfixCtrl(0),
158  _wcshdrCtrl(0),
159  _nReject(0),
160  _coordSystem(static_cast<afwCoord::CoordSystem>(-1))
161 {
163  initWcsLib(crval, crpix, CD,
164  ctype1, ctype2,
165  equinox, raDecSys,
166  cunits1, cunits2);
167  _initWcs();
168 }
169 
170 
177  class HeaderAccess {
178  public:
180  CONST_PTR(lsst::daf::base::PropertySet) const& toRead() { return _constHeader; }
182  PTR(lsst::daf::base::PropertySet) const& toWrite() {
183  if (!_hackHeader) {
184  _hackHeader = _constHeader->deepCopy();
185  _constHeader = _hackHeader;
186  }
187  return _hackHeader;
188  }
189 
191  HeaderAccess(CONST_PTR(lsst::daf::base::PropertySet) const& header) :
192  _constHeader(header), _hackHeader() {}
193 
194  private:
196  PTR(lsst::daf::base::PropertySet) _hackHeader;
197  };
198 
199  HeaderAccess access(header);
200 
201  // Some headers (e.g. SDSS ones from FNAL) have EQUINOX as a string. Fix this,
202  // as wcslib 4.4.4 refuses to handle it
203  {
204  std::string const& key = "EQUINOX";
205  if (access.toRead()->exists(key) && access.toRead()->typeOf(key) == typeid(std::string)) {
206  double equinox = ::atof(access.toRead()->getAsString(key).c_str());
207  access.toWrite()->set(key, equinox);
208  }
209  }
210 
211  //Check header isn't empty
212  int nCards = lsst::afw::formatters::countFitsHeaderCards(access.toRead());
213  if (nCards <= 0) {
214  string msg = "Could not parse FITS WCS: no header cards found";
215  throw LSST_EXCEPT(except::InvalidParameterError, msg);
216  }
217 
218  // Scamp produces PVi_xx header cards that are inconsistent with WCS Paper 2
219  // and cause WCSLib to choke. Aggressively, rename all PV keywords to X_PV
220  for (int j=1; j<3; j++) {
221  for (int i=0; i<=99; i++) {
222  std::string key = (boost::format("PV%i_%i") % j % i).str();
223  if (!access.toRead()->exists(key)) {
224  break;
225  }
226  double val = access.toRead()->getAsDouble(key);
227  access.toWrite()->remove(key);
228  access.toWrite()->add("X_"+key, val);
229  }
230  }
231 
232  //While the standard does not insist on CRVAL and CRPIX being present, it
233  //is almost certain their absence indicates a problem.
234  //Check for CRPIX
235  if( !access.toRead()->exists("CRPIX1") && !access.toRead()->exists("CRPIX1a")) {
236  string msg = "Neither CRPIX1 not CRPIX1a found";
237  throw LSST_EXCEPT(except::InvalidParameterError, msg);
238  }
239 
240  if( !access.toRead()->exists("CRPIX2") && !access.toRead()->exists("CRPIX2a")) {
241  string msg = "Neither CRPIX2 not CRPIX2a found";
242  throw LSST_EXCEPT(except::InvalidParameterError, msg);
243  }
244 
245  //And the same for CRVAL
246  if( !access.toRead()->exists("CRVAL1") && !access.toRead()->exists("CRVAL1a")) {
247  string msg = "Neither CRVAL1 not CRVAL1a found";
248  throw LSST_EXCEPT(except::InvalidParameterError, msg);
249  }
250 
251  if( !access.toRead()->exists("CRVAL2") && !access.toRead()->exists("CRVAL2a")) {
252  string msg = "Neither CRVAL2 not CRVAL2a found";
253  throw LSST_EXCEPT(except::InvalidParameterError, msg);
254  }
255  /*
256  * According to Greisen and Calabretta (A&A 395, 1061–1075 (2002)) it's illegal to mix PCi_j and CDi_j
257  * headers; unfortunately Subaru puts both in its headers. It actually uses PC001002 instead of PC1_2
258  * (dating to a proposed FITS standard from 1996) and at least sometimes fails to include CDELT[12],
259  * so the CD and PC matrices are inconsistent
260  *
261  * If we detect any part of a CD matrix, delete all PC matrices
262  */
263  if(access.toRead()->exists("CD1_1") || access.toRead()->exists("CD1_2") ||
264  access.toRead()->exists("CD2_1") || access.toRead()->exists("CD2_2")) {
265  for (int i = 1; i <= 2; ++i) {
266  for (int j = 1; j <= 2; ++j) {
267  std::string key = (boost::format("PC%i_%i") % j % i).str();
268  if (access.toRead()->exists(key)) {
269  double const val = access.toRead()->getAsDouble(key);
270  access.toWrite()->remove(key);
271  access.toWrite()->add("X_" + key, val);
272  }
273 
274  key = (boost::format("PC%03d%03d") % j % i).str();
275  if (access.toRead()->exists(key)) {
276  double const val = access.toRead()->getAsDouble(key);
277  access.toWrite()->remove(key);
278  access.toWrite()->add("X_" + key, val);
279  }
280  }
281  }
282  }
283 
284  //Pass the header into wcslib's formatter to extract & setup the Wcs. First need
285  //to convert to a C style string, so the compile doesn't complain about constness
286  std::string metadataStr = lsst::afw::formatters::formatFitsProperties(access.toRead());
287  // We own the data, and wcslib is slack about constness, so no qualms with casting away const
288  char *hdrString = const_cast<char*>(metadataStr.c_str());
289  //printf("wcspih string:\n%s\n", hdrString);
290 
291  nCards = lsst::afw::formatters::countFitsHeaderCards(access.toRead()); // we may have dropped some
292  int pihStatus = wcspih(hdrString, nCards, _relax, _wcshdrCtrl, &_nReject, &_nWcsInfo, &_wcsInfo);
293 
294  if (pihStatus != 0) {
295  throw LSST_EXCEPT(except::RuntimeError,
296  (boost::format("Could not parse FITS WCS: wcspih status = %d (%s)") %
297  pihStatus % wcs_errmsg[pihStatus]).str());
298  }
299 
300  //Run wcsfix on _wcsInfo to try and fix any problems it knows about.
301  const int *naxes = NULL; // should be {NAXIS1, NAXIS2, ...} to check cylindrical projections
302  int stats[NWCSFIX]; // status returns from wcsfix
303  int fixStatus = wcsfix(_wcsfixCtrl, naxes, _wcsInfo, stats);
304  if (fixStatus != 0) {
305  std::stringstream errStream;
306  errStream << "Could not parse FITS WCS: wcsfix failed " << std::endl;
307  for (int ii = 0; ii < NWCSFIX; ++ii) {
308  if (stats[ii] >= 0) {
309  errStream << "\t" << ii << ": " << stats[ii] << " " << wcsfix_errmsg[stats[ii]] << std::endl;
310  } else {
311  errStream << "\t" << ii << ": " << stats[ii] << std::endl;
312  }
313  }
314  }
315 
316  //The Wcs standard requires a default value for RADESYS if the keyword
317  //doesn't exist in header, but wcslib doesn't set it. So we do so here. This code
318  //conforms to Calabretta & Greisen 2002 \S 3.1
319  if (!(access.toRead()->exists("RADESYS") || access.toRead()->exists("RADESYSa"))) {
320  // If RADECSYS exists, use that (counter to Calabretta & Greisen 2002 \S 3.1, but commonly used).
321  // If equinox exist and < 1984, use FK4. If >= 1984, use FK5
322  if (access.toRead()->exists("RADECSYS")) {
323  strncpy(_wcsInfo->radesys, access.toRead()->getAsString("RADECSYS").c_str(), STRLEN);
324  } else if (access.toRead()->exists("EQUINOX") || access.toRead()->exists("EQUINOXa")) {
325  std::string const EQUINOX = access.toRead()->exists("EQUINOX") ? "EQUINOX" : "EQUINOXa";
326  double const equinox = access.toRead()->getAsDouble(EQUINOX);
327  if(equinox < 1984) {
328  strncpy(_wcsInfo->radesys, "FK4", STRLEN);
329  } else {
330  strncpy(_wcsInfo->radesys, "FK5", STRLEN);
331  }
332  } else {
333  //If Equinox doesn't exist, default to ICRS
334  strncpy(_wcsInfo->radesys, "ICRS", STRLEN);
335  }
336  }
337  // strip trailing whitespace
338  {
339  for(int i = strlen(_wcsInfo->radesys) - 1; i >= 0; i--) {
340  if (isspace(_wcsInfo->radesys[i])) {
341  _wcsInfo->radesys[i] = '\0';
342  }
343  }
344  }
345  //
346  // If there are no CDi_j cards in the header, set CDi_j from PCi_j
347  // CDi_j == CDELTi*PCi_j
348  //
349  if ((_wcsInfo->altlin & 2) == 0) { // no CDi_j cards were found in the header
350  double const *cdelt = _wcsInfo->cdelt;
351  double const *pc = _wcsInfo->pc;
352  double *cd = _wcsInfo->cd;
353 
354  cd[0] = cdelt[0]*pc[0]; // 1_1
355  cd[1] = cdelt[0]*pc[1]; // 1_2
356  cd[2] = cdelt[1]*pc[2]; // 2_1
357  cd[3] = cdelt[1]*pc[3]; // 2_2
358  }
359 }
360 
371 void Wcs::initWcsLib(GeomPoint const & crval, GeomPoint const & crpix, Eigen::Matrix2d const & CD,
372  std::string const & ctype1, std::string const & ctype2,
373  double equinox, std::string const & raDecSys,
374  std::string const & cunits1, std::string const & cunits2) {
375 
376  //Check CD is a valid size
377  if( (CD.rows() != 2) || (CD.cols() != 2) ) {
378  string msg = "CD is not a 2x2 matrix";
379  throw LSST_EXCEPT(except::InvalidParameterError, msg);
380  }
381 
382  //Check that cunits are legitimate values
383  bool isValid = (cunits1 == "deg");
384  isValid |= (cunits1 == "arcmin");
385  isValid |= (cunits1 == "arcsec");
386  isValid |= (cunits1 == "mas");
387 
388  if (!isValid) {
389  string msg = "CUNITS1 must be one of {deg|arcmin|arcsec|mas}";
390  throw LSST_EXCEPT(except::InvalidParameterError, msg);
391  }
392 
393  isValid = (cunits2 == "deg");
394  isValid |= (cunits2 == "arcmin");
395  isValid |= (cunits2 == "arcsec");
396  isValid |= (cunits2 == "mas");
397 
398  if (!isValid) {
399  string msg = "CUNITS2 must be one of {deg|arcmin|arcsec|mas}";
400  throw LSST_EXCEPT(except::InvalidParameterError, msg);
401  }
402 
403  //Initialise the wcs struct
404  _wcsInfo = static_cast<struct wcsprm *>(malloc(sizeof(struct wcsprm)));
405  if (_wcsInfo == NULL) {
406  throw LSST_EXCEPT(except::MemoryError, "Cannot allocate WCS info");
407  }
408 
409  _wcsInfo->flag = -1;
410  int status = wcsini(true, 2, _wcsInfo); //2 indicates a naxis==2, a two dimensional image
411  if(status != 0) {
412  throw LSST_EXCEPT(except::MemoryError,
413  (boost::format("Failed to allocate memory with wcsini. Status %d: %s") %
414  status % wcs_errmsg[status] ).str());
415  }
416 
417  //Set crval, crpix and CD. Internally to the class, we use fits units for consistency with
418  //wcslib.
419  _wcsInfo->crval[0] = crval.getX();
420  _wcsInfo->crval[1] = crval.getY();
421  _wcsInfo->crpix[0] = crpix.getX() + lsstToFitsPixels;
422  _wcsInfo->crpix[1] = crpix.getY() + lsstToFitsPixels;
423 
424  //Set the CD matrix
425  for (int i=0; i<2; ++i) {
426  for (int j=0; j<2; ++j) {
427  _wcsInfo->cd[(2*i) + j] = CD(i,j);
428  }
429  }
430 
431  //Specify that we have a CD matrix, but no PC or CROTA
432  _wcsInfo->altlin = 2;
433  _wcsInfo->flag = 0; //values have been updated
434 
435  //This is a work around for what I think is a bug in wcslib. ->types is neither
436  //initialised or set to NULL by default, so if I try to delete a Wcs object,
437  //wcslib then attempts to free non-existent space, and the code can crash.
438  _wcsInfo->types = NULL;
439 
440  //Set the coordinate system
441  strncpy(_wcsInfo->ctype[0], ctype1.c_str(), STRLEN);
442  strncpy(_wcsInfo->ctype[1], ctype2.c_str(), STRLEN);
443  strncpy(_wcsInfo->radesys, raDecSys.c_str(), STRLEN);
444  _wcsInfo->equinox = equinox;
445 
446  //Set the units
447  strncpy(_wcsInfo->cunit[0], cunits1.c_str(), STRLEN);
448  strncpy(_wcsInfo->cunit[1], cunits2.c_str(), STRLEN);
449 
450  _nWcsInfo = 1; //Specify that we have only one coordinate representation
451 
452  //Tell wcslib that we are need to set up internal values
453  status=wcsset(_wcsInfo);
454  if(status != 0) {
455  throw LSST_EXCEPT(except::RuntimeError,
456  (boost::format("Failed to setup wcs structure with wcsset. Status %d: %s") %
457  status % wcs_errmsg[status] ).str());
458 
459  }
460 }
461 
462 
464 Wcs::Wcs(afwImg::Wcs const & rhs) :
465  daf::base::Citizen(typeid(this)),
466  _wcsInfo(NULL),
467  _nWcsInfo(rhs._nWcsInfo),
468  _relax(rhs._relax),
469  _wcsfixCtrl(rhs._wcsfixCtrl),
470  _wcshdrCtrl(rhs._wcshdrCtrl),
471  _nReject(rhs._nReject),
472  _coordSystem(static_cast<afwCoord::CoordSystem>(-1))
473 {
474 
475  if (rhs._nWcsInfo > 0) {
476  _wcsInfo = static_cast<struct wcsprm *>(calloc(rhs._nWcsInfo, sizeof(struct wcsprm)));
477  if (_wcsInfo == NULL) {
478  throw LSST_EXCEPT(lsst::pex::exceptions::MemoryError, "Cannot allocate WCS info");
479  }
480 
481  _wcsInfo->flag = -1;
482  int alloc = 1; //Unconditionally allocate memory when calling
483  for (int i = 0; i != rhs._nWcsInfo; ++i) {
484  int status = wcscopy(alloc, &rhs._wcsInfo[i], &_wcsInfo[i]);
485  if (status != 0) {
486  wcsvfree(&i, &_wcsInfo);
487  throw LSST_EXCEPT(lsst::pex::exceptions::MemoryError,
488  (boost::format("Could not copy WCS: wcscopy status = %d : %s") %
489  status % wcs_errmsg[status]).str());
490  }
491  }
492  }
493  _initWcs();
494 }
495 
496 bool Wcs::operator==(Wcs const & other) const {
497  if (&other == this) return true;
498  // We do a bidirectional test with a virtual member function in case one of us is a derived
499  // class with members we don't know about here.
500  // This is not the most efficient possible implementation, but I think it's the easiest one
501  // with which to ensure correctness, and I think that's more important in this case.
502  return this->_isSubset(other) && other._isSubset(*this);
503 }
504 
505 // convenience functions and a macro for implementing isSubset
506 namespace {
507 
508 inline bool compareArrays(double const * a, double const * b, int n) {
509  for (int i = 0; i < n; ++i) if (a[i] != b[i]) return false;
510  return true;
511 }
512 
513 template <typename T>
514 inline bool compareStringArrays(T a, T b, int n) {
515  for (int i = 0; i < n; ++i) if (strcmp(a[i], b[i]) != 0) return false;
516  return true;
517 }
518 
519 #define CHECK_NULLS(a, b) \
520  do { \
521  if ((a) == NULL) { \
522  if ((b) == NULL) return true; \
523  return false; \
524  } \
525  if ((b) == NULL) return false; \
526  } while (false)
527 
528 } // anonymous
529 
530 bool Wcs::_isSubset(Wcs const & rhs) const {
532  CHECK_NULLS(_wcsInfo->crval, rhs._wcsInfo->crval);
533  CHECK_NULLS(_wcsInfo->cd, rhs._wcsInfo->cd);
534  CHECK_NULLS(_wcsInfo->crpix, rhs._wcsInfo->crpix);
535  CHECK_NULLS(_wcsInfo->cunit, rhs._wcsInfo->cunit);
536  CHECK_NULLS(_wcsInfo->ctype, rhs._wcsInfo->ctype);
537  return _nWcsInfo == rhs._nWcsInfo &&
538  _coordSystem == rhs._coordSystem &&
539  _wcsInfo->naxis == rhs._wcsInfo->naxis &&
540  _wcsInfo->equinox == rhs._wcsInfo->equinox &&
541  _wcsInfo->altlin == rhs._wcsInfo->altlin &&
542  compareArrays(_wcsInfo->crval, rhs._wcsInfo->crval, 2) &&
543  compareArrays(_wcsInfo->crpix, rhs._wcsInfo->crpix, 2) &&
544  compareArrays(_wcsInfo->cd, rhs._wcsInfo->cd, 4) &&
545  compareStringArrays(_wcsInfo->cunit, rhs._wcsInfo->cunit, 2) &&
546  compareStringArrays(_wcsInfo->ctype, rhs._wcsInfo->ctype, 2) &&
548  _wcsInfo->crval[1] * afwGeom::degrees) ==
549  rhs.skyToPixel(_wcsInfo->crval[0] * afwGeom::degrees,
550  _wcsInfo->crval[1] * afwGeom::degrees) &&
551  *pixelToSky(_wcsInfo->crpix[0], _wcsInfo->crpix[1]) ==
552  *rhs.pixelToSky(_wcsInfo->crpix[0], _wcsInfo->crpix[1]);
553 }
554 
556  if (_wcsInfo != NULL) {
557  wcsvfree(&_nWcsInfo, &_wcsInfo);
558  }
559 }
560 
561 
562 Wcs::Ptr Wcs::clone(void) const {
563  return Wcs::Ptr(new Wcs(*this));
564 }
565 
566 //
567 // Accessors
568 //
569 
571  assert(_wcsInfo);
572  return makeCorrectCoord(_wcsInfo->crval[0] * afwGeom::degrees, _wcsInfo->crval[1] * afwGeom::degrees);
573 }
574 
576  assert(_wcsInfo);
577  //Convert from fits units back to lsst units
578  double p1 = _wcsInfo->crpix[0] + fitsToLsstPixels;
579  double p2 = _wcsInfo->crpix[1] + fitsToLsstPixels;
580  return afwGeom::Point2D(p1, p2);
581 }
582 
583 Eigen::Matrix2d Wcs::getCDMatrix() const {
584  assert(_wcsInfo);
585  int const naxis = _wcsInfo->naxis;
586 
587  //If naxis != 2, I'm not sure if any of what follows is correct
588  assert(naxis == 2);
589 
590  Eigen::Matrix2d C;
591 
592  for (int i=0; i< naxis; ++i){
593  for (int j=0; j<naxis; ++j) {
594  C(i,j) = _wcsInfo->cd[ (i*naxis) + j ];
595  }
596  }
597 
598  return C;
599 }
600 
601 void Wcs::flipImage(int flipLR, int flipTB, afwGeom::Extent2I dimensions) const {
602  assert(_wcsInfo);
603 
604  int const naxis = _wcsInfo->naxis;
605 
606  //If naxis != 2, I'm not sure if any of what follows is correct
607  assert(naxis == 2);
608  if (flipLR) {
609  _wcsInfo->cd[0] = -_wcsInfo->cd[0];
610  _wcsInfo->cd[2] = -_wcsInfo->cd[2];
611  _wcsInfo->crpix[0] = -_wcsInfo->crpix[0] + dimensions.getX();
612  }
613  if (flipTB) {
614  _wcsInfo->cd[1] = -_wcsInfo->cd[1];
615  _wcsInfo->cd[3] = -_wcsInfo->cd[3];
616  _wcsInfo->crpix[1] = -_wcsInfo->crpix[1]+dimensions.getY();
617  }
618 
619  // tells libwcs to invalidate cached data, since transformation has been modified
620  _wcsInfo->flag = 0;
621 }
622 
624  assert(_wcsInfo);
625 
626  while (nQuarter < 0 ) {
627  nQuarter += 4;
628  }
629 
630 
631  int const naxis = _wcsInfo->naxis;
632 
633  //If naxis != 2, I'm not sure if any of what follows is correct
634  assert(naxis == 2);
635  double a = _wcsInfo->cd[0];
636  double b = _wcsInfo->cd[1];
637  double c = _wcsInfo->cd[2];
638  double d = _wcsInfo->cd[3];
639  double crpx = _wcsInfo->crpix[0];
640  double crpy = _wcsInfo->crpix[1];
641  switch (nQuarter%4) {
642  case 0:
643  break;
644  case 1:
645  _wcsInfo->cd[0] = -b;
646  _wcsInfo->cd[1] = a;
647  _wcsInfo->cd[2] = -d;
648  _wcsInfo->cd[3] = c;
649  _wcsInfo->crpix[0] = -crpy + dimensions.getY();
650  _wcsInfo->crpix[1] = crpx;
651  break;
652  case 2:
653  _wcsInfo->cd[0] = -a;
654  _wcsInfo->cd[1] = -b;
655  _wcsInfo->cd[2] = -c;
656  _wcsInfo->cd[3] = -d;
657  _wcsInfo->crpix[0] = -crpx + dimensions.getX();
658  _wcsInfo->crpix[1] = -crpy + dimensions.getY();
659  break;
660  case 3:
661  _wcsInfo->cd[0] = b;
662  _wcsInfo->cd[1] = -a;
663  _wcsInfo->cd[2] = d;
664  _wcsInfo->cd[3] = -c;
665  _wcsInfo->crpix[0] = crpy;
666  _wcsInfo->crpix[1] = -crpx + dimensions.getX();
667  break;
668  }
669 
670  // tells libwcs to invalidate cached data, since transformation has been modified
671  _wcsInfo->flag = 0;
672 }
673 
676 }
677 
678 bool Wcs::isFlipped() const {
679  // We calculate the determinant of the CD (i.e the rotation and scaling)
680  // matrix. If this determinant is positive, then the image can be rotated
681  // to a position where increasing the right ascension and declination
682  // increases the horizontal and vertical pixel position. In this case the
683  // image is flipped.
684  assert(_wcsInfo);
685  double det = (_wcsInfo->cd[0] * _wcsInfo->cd[3]) - (_wcsInfo->cd[1] * _wcsInfo->cd[2]);
686 
687  if (det == 0) {
688  throw(LSST_EXCEPT(except::RuntimeError, "Wcs CD matrix is singular"));
689  }
690 
691  return (det > 0);
692 }
693 
694 static double square(double x) {
695  return x*x;
696 }
697 
698 double Wcs::pixArea(GeomPoint pix00
699  ) const {
700  //
701  // Figure out the (0, 0), (0, 1), and (1, 0) ra/dec coordinates of the corners of a square drawn in pixel
702  // It'd be better to centre the square at sky00, but that would involve another conversion between sky and
703  // pixel coordinates so I didn't bother
704  //
705  const double side = 1; // length of the square's sides in pixels
706 
707  // Work in 3-space to avoid RA wrapping and pole issues.
708  afwGeom::Point3D v0 = pixelToSky(pix00)->getVector();
709 
710  // Step by "side" in x and y pixel directions...
711  GeomPoint px(pix00);
712  GeomPoint py(pix00);
713  px.shift(afwGeom::Extent2D(side, 0));
714  py.shift(afwGeom::Extent2D(0, side));
715  // Push the points through the WCS, and find difference in 3-space.
716  afwGeom::Extent3D dx = pixelToSky(px)->getVector() - v0;
717  afwGeom::Extent3D dy = pixelToSky(py)->getVector() - v0;
718 
719  // Compute |cross product| = area of parallelogram with sides dx,dy
720  // FIXME -- this is slightly incorrect -- it's making the small-angle
721  // approximation, taking the distance *through* the unit sphere
722  // rather than over its surface.
723  // This is in units of ~radians^2
724  double area = sqrt(square(dx[1]*dy[2] - dx[2]*dy[1]) +
725  square(dx[2]*dy[0] - dx[0]*dy[2]) +
726  square(dx[0]*dy[1] - dx[1]*dy[0]));
727 
728  return area / square(side) * square(180. / afwGeom::PI);
729 }
730 
732  return sqrt(pixArea(getPixelOrigin())) * afwGeom::degrees;
733 }
734 
735 /*
736  * Worker routine for skyToPixel
737  */
738 GeomPoint Wcs::skyToPixelImpl(afwGeom::Angle sky1, // RA (or, more generally, longitude)
739  afwGeom::Angle sky2 // Dec (or latitude)
740  ) const {
741  assert(_wcsInfo);
742 
743  double skyTmp[2];
744  double imgcrd[2];
745  double phi, theta;
746  double pixTmp[2];
747  /*
748  printf("_skyCoordsReversed: %c\n", (_skyCoordsReversed ? 'T' : 'F'));
749  printf("wcsinfo.lat: %i, lng: %i\n", _wcsInfo->lat, _wcsInfo->lng);
750  */
751  // WCSLib is smart enough to notice and handle crazy SDSS CTYPE1 = DEC--TAN,
752  // by recording the indices of the long and lat coordinates.
753  skyTmp[_wcsInfo->lng] = sky1.asDegrees();
754  skyTmp[_wcsInfo->lat] = sky2.asDegrees();
755 
756  int stat[1];
757  int status = 0;
758  status = wcss2p(_wcsInfo, 1, 2, skyTmp, &phi, &theta, imgcrd, pixTmp, stat);
759  if (status == 9) {
760  throw LSST_EXCEPT(except::DomainError,
761  (boost::format("sky coordinates %s, %s degrees is not valid for this WCS")
762  % sky1.asDegrees() % sky2.asDegrees()
763  ).str()
764  );
765  }
766  if (status > 0) {
767  throw LSST_EXCEPT(except::RuntimeError,
768  (boost::format("Error: wcslib returned a status code of %d at sky %s, %s deg: %s") %
769  status % sky1.asDegrees() % sky2.asDegrees() % wcs_errmsg[status]).str());
770  }
771 
772  // wcslib assumes 1-indexed coords
775 }
776 
779  return skyToPixelImpl(sky->getLongitude(), sky->getLatitude());
780 }
781 
782 
785  return coord.convert(_coordSystem);
786 }
787 
789  return skyToPixelImpl(sky1, sky2);
790 }
791 
793  assert(_wcsInfo);
794 
796  double skyTmp[2];
797  double imgcrd[2];
798  double phi, theta;
799  double pixTmp[2];
800 
801  /*
802  printf("skyToIWC: _coordSystem = %i\n", (int)_coordSystem);
803  printf("coord (%.3f, %.3f)\n", coord->getLongitude().asDegrees(), coord->getLatitude().asDegrees());
804  printf("->sky (%.3f, %.3f)\n", sky->getLongitude().asDegrees(), sky->getLatitude().asDegrees());
805  */
806 
807  skyTmp[_wcsInfo->lng] = sky->getLongitude().asDegrees();
808  skyTmp[_wcsInfo->lat] = sky->getLatitude() .asDegrees();
809 
810  //Estimate pixel coordinates
811  int stat[1];
812  int status = 0;
813  imgcrd[0] = imgcrd[1] = -1e6;
814  /*
815  printf(" skyTmp[] = (%.3f, %.3f)\n", skyTmp[0], skyTmp[1]);
816  printf(" _wcsInfo->lng,lat = %i, %i\n", _wcsInfo->lng, _wcsInfo->lat);
817  */
818  status = wcss2p(_wcsInfo, 1, 2, skyTmp, &phi, &theta, imgcrd, pixTmp, stat);
819  if (status > 0) {
820  throw LSST_EXCEPT(except::RuntimeError,
821  (boost::format("Error: wcslib returned a status code of %d at sky %s, %s deg: %s") %
822  status % skyTmp[0] % skyTmp[1] % wcs_errmsg[status]).str());
823  }
824  /*
825  printf("->iwc (%.3f, %.3f)\n", imgcrd[0], imgcrd[1]);
826  printf("-> pix (%.2f, %.2f)\n", pixTmp[0], pixTmp[1]);
827  afwCoord::Coord::Ptr crval = getSkyOrigin();
828  printf("(crval is (%.3f, %.3f))\n", crval->getLongitude().asDegrees(), crval->getLatitude().asDegrees());
829  */
830  return GeomPoint(imgcrd[0], imgcrd[1]);
831 }
832 
833 /*
834  * Worker routine for pixelToSky
835  */
836 void
837 Wcs::pixelToSkyImpl(double pixel1, double pixel2, afwGeom::Angle skyTmp[2]) const
838 {
839  assert(_wcsInfo);
840 
841  // wcslib assumes 1-indexed coordinates
842  double pixTmp[2] = { pixel1 - lsst::afw::image::PixelZeroPos + lsstToFitsPixels,
843  pixel2 - lsst::afw::image::PixelZeroPos + lsstToFitsPixels};
844  double imgcrd[2];
845  double phi, theta;
846 
847  double sky[2];
848  int status = 0;
849  status = wcsp2s(_wcsInfo, 1, 2, pixTmp, imgcrd, &phi, &theta, sky, &status);
850  if (status > 0) {
851  throw LSST_EXCEPT(except::RuntimeError,
852  (boost::format("Error: wcslib returned a status code of %d at pixel %s, %s: %s") %
853  status % pixel1 % pixel2 % wcs_errmsg[status]).str());
854  }
855  // FIXME -- _wcsInfo.lat, _wcsInfo.lng ?
856  skyTmp[0] = sky[0] * afwGeom::degrees;
857  skyTmp[1] = sky[1] * afwGeom::degrees;
858 }
859 
861  return pixelToSky(pixel.getX(), pixel.getY());
862 }
863 
864 CoordPtr Wcs::pixelToSky(double pixel1, double pixel2) const {
865  assert(_wcsInfo);
866 
867  afwGeom::Angle skyTmp[2];
868  pixelToSkyImpl(pixel1, pixel2, skyTmp);
869  return makeCorrectCoord(skyTmp[0], skyTmp[1]);
870 }
871 
872 void Wcs::pixelToSky(double pixel1, double pixel2, afwGeom::Angle& sky1, afwGeom::Angle& sky2) const {
873  afwGeom::Angle skyTmp[2];
874  // HACK -- we shouldn't need to initialize these -- pixelToSkyImpl() sets them unless an
875  // exception is thrown -- but be safe.
876  skyTmp[0] = 0. * afwGeom::radians;
877  skyTmp[1] = 0. * afwGeom::radians;
878  pixelToSkyImpl(pixel1, pixel2, skyTmp);
879  sky1 = skyTmp[0];
880  sky2 = skyTmp[1];
881 }
882 
886 
887  //Construct a coord object of the correct type
888  int const ncompare = 4; // we only care about type's first 4 chars
889  char *type = _wcsInfo->ctype[0];
890  char *radesys = _wcsInfo->radesys;
891  double equinox = _wcsInfo->equinox;
892 
893  if (strncmp(type, "RA--", ncompare) == 0) { // Our default. If it's often something else, consider
894  ; // using an tr1::unordered_map
895  if(strcmp(radesys, "ICRS") == 0) {
896  return afwCoord::makeCoord(afwCoord::ICRS, sky0, sky1);
897  }
898  if(strcmp(radesys, "FK5") == 0) {
899  return afwCoord::makeCoord(afwCoord::FK5, sky0, sky1, equinox);
900  } else {
901  throw LSST_EXCEPT(except::RuntimeError,
902  (boost::format("Can't create Coord object: Unrecognised radesys %s") %
903  radesys).str());
904  }
905 
906  } else if (strncmp(type, "GLON", ncompare) == 0) {
907  return afwCoord::makeCoord(afwCoord::GALACTIC, sky0, sky1);
908  } else if (strncmp(type, "ELON", ncompare) == 0) {
909  return afwCoord::makeCoord(afwCoord::ECLIPTIC, sky0, sky1, equinox);
910  } else if (strncmp(type, "DEC-", ncompare) == 0) {
911  //check for the case where the ctypes are swapped. Note how sky0 and sky1 are swapped as well
912 
913  //Our default
914  if(strcmp(radesys, "ICRS") == 0) {
915  return afwCoord::makeCoord(afwCoord::ICRS, sky1, sky0);
916  }
917  if(strcmp(radesys, "FK5") == 0) {
918  return afwCoord::makeCoord(afwCoord::FK5, sky1, sky0, equinox);
919  } else {
920  throw LSST_EXCEPT(except::RuntimeError,
921  (boost::format("Can't create Coord object: Unrecognised radesys %s") %
922  radesys).str());
923  }
924  } else if (strncmp(type, "GLAT", ncompare) == 0) {
925  return afwCoord::makeCoord(afwCoord::GALACTIC, sky1, sky0);
926  } else if (strncmp(type, "ELAT", ncompare) == 0) {
927  return afwCoord::makeCoord(afwCoord::ECLIPTIC, sky1, sky0, equinox);
928  } else {
929  //Give up in disgust
930  throw LSST_EXCEPT(except::RuntimeError,
931  (boost::format("Can't create Coord object: Unrecognised sys %s") %
932  type).str());
933  }
934 
935  //Can't get here
936  assert(0);
937 }
938 
939 
941  lsst::afw::coord::Coord const & coord,
943 ) const {
944  return linearizePixelToSkyInternal(skyToPixel(coord), coord, skyUnit);
945 }
947  GeomPoint const & pix,
949 ) const {
950  return linearizePixelToSkyInternal(pix, *pixelToSky(pix), skyUnit);
951 }
952 
953 /*
954  * Implementation for the overloaded public linearizePixelToSky methods, requiring both a pixel coordinate
955  * and the corresponding sky coordinate.
956  */
958  GeomPoint const & pix00,
959  lsst::afw::coord::Coord const & coord,
961 ) const {
962  //
963  // Figure out the (0, 0), (0, 1), and (1, 0) ra/dec coordinates of the corners of a square drawn in pixel
964  // It'd be better to centre the square at sky00, but that would involve another conversion between sky and
965  // pixel coordinates so I didn't bother
966  //
967  const double side = 10; // length of the square's sides in pixels
968  GeomPoint const sky00 = coord.getPosition(skyUnit);
969  typedef std::pair<lsst::afw::geom::Angle, lsst::afw::geom::Angle> AngleAngle;
970  AngleAngle const dsky10 = coord.getTangentPlaneOffset(*pixelToSky(pix00 + afwGeom::Extent2D(side, 0)));
971  AngleAngle const dsky01 = coord.getTangentPlaneOffset(*pixelToSky(pix00 + afwGeom::Extent2D(0, side)));
972 
973  Eigen::Matrix2d m;
974  m(0, 0) = dsky10.first.asAngularUnits(skyUnit)/side;
975  m(0, 1) = dsky01.first.asAngularUnits(skyUnit)/side;
976  m(1, 0) = dsky10.second.asAngularUnits(skyUnit)/side;
977  m(1, 1) = dsky01.second.asAngularUnits(skyUnit)/side;
978 
979  Eigen::Vector2d sky00v;
980  sky00v << sky00.getX(), sky00.getY();
981  Eigen::Vector2d pix00v;
982  pix00v << pix00.getX(), pix00.getY();
983  //return lsst::afw::geom::AffineTransform(m, lsst::afw::geom::Extent2D(sky00v - m * pix00v));
984  return lsst::afw::geom::AffineTransform(m, (sky00v - m * pix00v));
985 }
986 
988  lsst::afw::coord::Coord const & coord,
990 ) const {
991  return linearizeSkyToPixelInternal(skyToPixel(coord), coord, skyUnit);
992 }
993 
995  GeomPoint const & pix,
997 ) const {
998  return linearizeSkyToPixelInternal(pix, *pixelToSky(pix), skyUnit);
999 }
1000 
1006  GeomPoint const & pix00,
1007  lsst::afw::coord::Coord const & coord,
1009 ) const {
1010  lsst::afw::geom::AffineTransform inverse = linearizePixelToSkyInternal(pix00, coord, skyUnit);
1011  return inverse.invert();
1012 }
1013 
1015 {
1017 }
1018 
1019 // -------- table-based persistence -------------------------------------------------------------------------
1020 
1021 namespace lsst { namespace afw { namespace image {
1022 
1024 public:
1025 
1026  explicit WcsFactory(std::string const & name) : table::io::PersistableFactory(name) {}
1027 
1028  virtual PTR(table::io::Persistable) read(
1029  InputArchive const & archive,
1030  CatalogVector const & catalogs
1031  ) const;
1032 
1033 };
1034 
1035 namespace {
1036 
1037 // Read-only singleton struct containing the schema and keys that a simple Wcs is mapped
1038 // to in record persistence.
1039 struct WcsPersistenceHelper : private boost::noncopyable {
1050 
1051  static WcsPersistenceHelper const & get() {
1052  static WcsPersistenceHelper instance;
1053  return instance;
1054  };
1055 
1056 private:
1057  WcsPersistenceHelper() :
1058  schema(),
1059  crval(schema.addField< table::Point<double> >("crval", "celestial reference point")),
1060  crpix(schema.addField< table::Point<double> >("crpix", "pixel reference point")),
1061  cd(schema.addField< table::Array<double> >(
1062  "cd", "linear transform matrix, ordered (1_1, 2_1, 1_2, 2_2)", 4)),
1063  ctype1(schema.addField< std::string >("ctype1", "coordinate type", 72)),
1064  ctype2(schema.addField< std::string >("ctype2", "coordinate type", 72)),
1065  equinox(schema.addField< double >("equinox", "equinox of coordinates")),
1066  radesys(schema.addField< std::string >("radesys", "coordinate system for equinox", 72)),
1067  cunit1(schema.addField< std::string >("cunit1", "coordinate units", 72)),
1068  cunit2(schema.addField< std::string >("cunit2", "coordinate units", 72))
1069  {
1070  schema.getCitizen().markPersistent();
1071  }
1072 };
1073 
1074 std::string getWcsPersistenceName() { return "Wcs"; }
1075 
1076 WcsFactory registration(getWcsPersistenceName());
1077 
1078 } // anonymous
1079 
1080 std::string Wcs::getPersistenceName() const { return getWcsPersistenceName(); }
1081 
1082 std::string Wcs::getPythonModule() const { return "lsst.afw.image"; }
1083 
1084 void Wcs::write(OutputArchiveHandle & handle) const {
1085  WcsPersistenceHelper const & keys = WcsPersistenceHelper::get();
1086  afw::table::BaseCatalog catalog = handle.makeCatalog(keys.schema);
1087  PTR(afw::table::BaseRecord) record = catalog.addNew();
1088  record->set(keys.crval, getSkyOrigin()->getPosition(afw::geom::degrees));
1089  record->set(keys.crpix, getPixelOrigin());
1090  Eigen::Matrix2d cdIn = getCDMatrix();
1091  Eigen::Map<Eigen::Matrix2d> cdOut((*record)[keys.cd].getData());
1092  cdOut = cdIn;
1093  record->set(keys.ctype1, std::string(_wcsInfo[0].ctype[0]));
1094  record->set(keys.ctype2, std::string(_wcsInfo[0].ctype[1]));
1095  record->set(keys.equinox, _wcsInfo[0].equinox);
1096  record->set(keys.radesys, std::string(_wcsInfo[0].radesys));
1097  record->set(keys.cunit1, std::string(_wcsInfo[0].cunit[0]));
1098  record->set(keys.cunit2, std::string(_wcsInfo[0].cunit[1]));
1099  handle.saveCatalog(catalog);
1100 }
1101 
1102 bool Wcs::isPersistable() const {
1103  if (_wcsInfo[0].naxis != 2) return false;
1104  if (std::strcmp(_wcsInfo[0].cunit[0], "deg") != 0) return false;
1105  if (std::strcmp(_wcsInfo[0].cunit[1], "deg") != 0) return false;
1106  return true;
1107 }
1108 
1110  daf::base::Citizen(typeid(this)),
1111  _wcsInfo(NULL),
1112  _nWcsInfo(0),
1113  _relax(0),
1114  _wcsfixCtrl(0),
1115  _wcshdrCtrl(0),
1116  _nReject(0),
1117  _coordSystem(static_cast<afw::coord::CoordSystem>(-1))
1118 {
1119  WcsPersistenceHelper const & keys = WcsPersistenceHelper::get();
1120  if (!record.getSchema().contains(keys.schema)) {
1121  throw LSST_EXCEPT(
1122  afw::table::io::MalformedArchiveError,
1123  "Incorrect schema for Wcs persistence"
1124  );
1125  }
1126  _setWcslibParams();
1127  Eigen::Matrix2d cd = Eigen::Map<Eigen::Matrix2d const>(record[keys.cd].getData());
1128  initWcsLib(
1129  record.get(keys.crval), record.get(keys.crpix), cd,
1130  record.get(keys.ctype1), record.get(keys.ctype2),
1131  record.get(keys.equinox), record.get(keys.radesys),
1132  record.get(keys.cunit1), record.get(keys.cunit2)
1133  );
1134  _initWcs();
1135 }
1136 
1138 WcsFactory::read(InputArchive const & inputs, CatalogVector const & catalogs) const {
1139  WcsPersistenceHelper const & keys = WcsPersistenceHelper::get();
1140  LSST_ARCHIVE_ASSERT(catalogs.size() >= 1u);
1141  LSST_ARCHIVE_ASSERT(catalogs.front().size() == 1u);
1142  LSST_ARCHIVE_ASSERT(catalogs.front().getSchema() == keys.schema);
1143  PTR(Wcs) result(new Wcs(catalogs.front().front()));
1144  return result;
1145 }
1146 
1147 }}} // namespace lsst::afw::image
1148 
1149 // ----------------------------------------------------------------------------------------------------------
1150 
1151 //
1152 //Mutators
1153 //
1154 
1155 
1161 void Wcs::shiftReferencePixel(double dx, double dy) {
1162  assert(_wcsInfo);
1163  _wcsInfo->crpix[0] += dx;
1164  _wcsInfo->crpix[1] += dy;
1165 
1166  // tells libwcs to invalidate cached data, since transformation has been modified
1167  _wcsInfo->flag = 0;
1168 }
1169 
1170 /************************************************************************************************************/
1171 /*
1172  * Now WCSA, pixel coordinates, but allowing for X0 and Y0
1173  */
1174 namespace lsst {
1175 namespace afw {
1176 namespace image {
1177 
1178 namespace detail {
1179 
1184 createTrivialWcsAsPropertySet(std::string const& wcsName,
1185  int const x0,
1186  int const y0
1187  ) {
1189 
1190  wcsMetaData->set("CRVAL1" + wcsName, x0, "Column pixel of Reference Pixel");
1191  wcsMetaData->set("CRVAL2" + wcsName, y0, "Row pixel of Reference Pixel");
1192  wcsMetaData->set("CRPIX1" + wcsName, 1, "Column Pixel Coordinate of Reference");
1193  wcsMetaData->set("CRPIX2" + wcsName, 1, "Row Pixel Coordinate of Reference");
1194  wcsMetaData->set("CTYPE1" + wcsName, "LINEAR", "Type of projection");
1195  wcsMetaData->set("CTYPE2" + wcsName, "LINEAR", "Type of projection");
1196  wcsMetaData->set("CUNIT1" + wcsName, "PIXEL", "Column unit");
1197  wcsMetaData->set("CUNIT2" + wcsName, "PIXEL", "Row unit");
1198 
1199  return wcsMetaData;
1200 }
1207 afwGeom::Point2I getImageXY0FromMetadata(std::string const& wcsName,
1208  lsst::daf::base::PropertySet *metadata
1209  ) {
1210 
1211  int x0 = 0; // Our value of X0
1212  int y0 = 0; // Our value of Y0
1213 
1214  try {
1215  //
1216  // Only use WCS if CRPIX[12] == 1 and CRVAL[12] is present
1217  //
1218  if (metadata->getAsDouble("CRPIX1" + wcsName) == 1 &&
1219  metadata->getAsDouble("CRPIX2" + wcsName) == 1) {
1220 
1221  x0 = metadata->getAsInt("CRVAL1" + wcsName);
1222  y0 = metadata->getAsInt("CRVAL2" + wcsName);
1223  //
1224  // OK, we've got it. Remove it from the header
1225  //
1226  metadata->remove("CRVAL1" + wcsName);
1227  metadata->remove("CRVAL2" + wcsName);
1228  metadata->remove("CRPIX1" + wcsName);
1229  metadata->remove("CRPIX2" + wcsName);
1230  metadata->remove("CTYPE1" + wcsName);
1231  metadata->remove("CTYPE1" + wcsName);
1232  metadata->remove("CUNIT1" + wcsName);
1233  metadata->remove("CUNIT2" + wcsName);
1234  }
1235  } catch(lsst::pex::exceptions::NotFoundError &) {
1236  ; // OK, not present
1237  }
1238 
1239  return afwGeom::Point2I(x0, y0);
1240 }
1241 
1251  CONST_PTR(Wcs) const& wcs
1252  )
1253 {
1254  PTR(lsst::daf::base::PropertySet) wcsMetadata = wcs->getFitsMetadata();
1255  std::vector<std::string> paramNames = wcsMetadata->paramNames();
1256  paramNames.push_back("CDELT1");
1257  paramNames.push_back("CDELT2");
1258  paramNames.push_back("LTV1");
1259  paramNames.push_back("LTV2");
1260  paramNames.push_back("PC001001");
1261  paramNames.push_back("PC001002");
1262  paramNames.push_back("PC002001");
1263  paramNames.push_back("PC002002");
1264  for (std::vector<std::string>::const_iterator ptr = paramNames.begin(); ptr != paramNames.end(); ++ptr) {
1265  metadata->remove(*ptr);
1266  }
1267 
1268  return 0; // would be ncard if remove returned a status
1269 }
1270 
1271 
1272 }}}}
1273 
1274 
1275 
1276 // -------------------------------------------------------------------------------------------------
1277 //
1278 // XYTransformFromWcsPair
1279 
1280 
1282  : XYTransform(), _dst(dst), _src(src)
1283 { }
1284 
1285 
1287 {
1288  return boost::make_shared<XYTransformFromWcsPair>(_dst->clone(), _src->clone());
1289 }
1290 
1291 
1293 {
1294  //
1295  // TODO there is an alternate version of pixelToSky() which is designated for the
1296  // "knowledgeable user in need of performance". This is probably better, but first I need
1297  // to understand exactly which checks are needed (e.g. I think we need to check by hand
1298  // that both Wcs's use the same celestial coordinate system)
1299  //
1300  PTR(afw::coord::Coord) x = _src->pixelToSky(pixel);
1301  return _dst->skyToPixel(*x);
1302 }
1303 
1305 {
1306  PTR(afw::coord::Coord) x = _dst->pixelToSky(pixel);
1307  return _src->skyToPixel(*x);
1308 }
1309 
1311 {
1312  // just swap src, dst
1313  return boost::make_shared<XYTransformFromWcsPair> (_src, _dst);
1314 }
Defines the fields and offsets for a table.
Definition: Schema.h:46
int _wcsfixCtrl
Do potentially unsafe translations of non-standard unit strings? 0/1 = no/yes.
Definition: Wcs.h:393
geom::Point2I getImageXY0FromMetadata(std::string const &wcsName, lsst::daf::base::PropertySet *metadata)
Definition: Wcs.cc:1207
const int lsstToFitsPixels
Definition: Wcs.cc:78
Interface for WcsFormatter class.
boost::shared_ptr< Coord > Ptr
Definition: Coord.h:72
#define CONST_PTR(...)
Definition: base.h:47
Schema getSchema() const
Return the Schema that holds this record&#39;s fields and keys.
Definition: BaseRecord.h:48
virtual std::string getPersistenceName() const
Return the unique name used to persist this object and look up its factory.
Definition: Wcs.cc:1080
double dx
Definition: ImageUtils.cc:90
Image utility functions.
boost::shared_ptr< Wcs const > _src
Definition: Wcs.h:446
BaseCatalog makeCatalog(Schema const &schema)
Return a new, empty catalog with the given schema.
geom::AffineTransform linearizePixelToSky(coord::Coord const &coord, geom::AngleUnit skyUnit=geom::degrees) const
Return the local linear approximation to Wcs::pixelToSky at a point given in sky coordinates.
Definition: Wcs.cc:940
WcsFactory(std::string const &name)
Definition: Wcs.cc:1026
table::Key< table::Array< double > > cd
Definition: Wcs.cc:1043
An object passed to Persistable::write to allow it to persist itself.
AngleUnit const radians
constant with units of radians
Definition: Angle.h:91
Formatting utilities.
AffineTransform const invert() const
virtual void pixelToSkyImpl(double pixel1, double pixel2, geom::Angle skyTmp[2]) const
Definition: Wcs.cc:837
coord::CoordSystem _coordSystem
Definition: Wcs.h:396
void _setWcslibParams()
Definition: Wcs.cc:68
A custom container class for records, based on std::vector.
Definition: Catalog.h:94
Class for storing ordered metadata with comments.
Definition: PropertyList.h:81
afw::table::Key< afw::table::Point< int > > dimensions
boost::shared_ptr< Wcs > Ptr
Definition: Wcs.h:113
A base class for factory classes used to reconstruct objects from records.
Definition: Persistable.h:231
geom::Point2D skyToPixel(geom::Angle sky1, geom::Angle sky2) const
Convert from sky coordinates (e.g. RA/dec) to pixel positions.
Definition: Wcs.cc:788
afw::coord::Coord::Ptr convertCoordToSky(coord::Coord const &coord) const
Definition: Wcs.cc:784
std::pair< lsst::afw::geom::Angle, lsst::afw::geom::Angle > getTangentPlaneOffset(Coord const &c) const
Get the offset on the tangent plane.
Definition: Coord.cc:747
virtual bool _isSubset(Wcs const &other) const
Definition: Wcs.cc:530
XYTransformFromWcsPair: An XYTransform obtained by putting two Wcs objects &quot;back to back&quot;...
Definition: Wcs.h:431
lsst::afw::geom::Point2D GeomPoint
Definition: Wcs.cc:60
afw::coord::Coord::Ptr makeCorrectCoord(geom::Angle sky0, geom::Angle sky1) const
Given a sky position, use the values stored in ctype and radesys to return the correct sub-class of C...
Definition: Wcs.cc:885
void initWcsLibFromFits(boost::shared_ptr< lsst::daf::base::PropertySet const > const &fitsMetadata)
Parse a fits header, extract the relevant metadata and create a Wcs object.
Definition: Wcs.cc:172
boost::shared_ptr< PropertyList > Ptr
Definition: PropertyList.h:84
void initWcsLib(geom::Point2D const &crval, geom::Point2D const &crpix, Eigen::Matrix2d const &CD, std::string const &ctype1, std::string const &ctype2, double equinox, std::string const &raDecSys, std::string const &cunits1, std::string const &cunits2)
Manually initialise a wcs struct using values passed by the constructor.
Definition: Wcs.cc:371
static lsst::daf::base::PropertyList::Ptr generatePropertySet(lsst::afw::image::Wcs const &wcs)
virtual geom::AffineTransform linearizeSkyToPixelInternal(geom::Point2D const &pix, coord::Coord const &coord, geom::AngleUnit skyUnit) const
Definition: Wcs.cc:1005
table::Key< std::string > radesys
Definition: Wcs.cc:1047
geom::Point2D skyToIntermediateWorldCoord(coord::Coord const &coord) const
Convert from sky coordinates (e.g. RA/dec) to intermediate world coordinates.
Definition: Wcs.cc:792
virtual ~Wcs()
Definition: Wcs.cc:555
virtual Ptr clone(void) const
Definition: Wcs.cc:562
Implementation of the WCS standard for a any projection.
Definition: Wcs.h:107
#define PTR(...)
Definition: base.h:41
Point< int, 2 > Point2I
Definition: PSF.h:39
double dy
Definition: ImageUtils.cc:90
boost::shared_ptr< coord::Coord > pixelToSky(double pix1, double pix2) const
Convert from pixel position to sky coordinates (e.g. RA/dec)
Definition: Wcs.cc:864
lsst::daf::base::PropertyList PropertyList
Definition: Wcs.cc:58
table::Key< std::string > ctype2
Definition: Wcs.cc:1045
A class used to convert scalar POD types such as double to Angle.
Definition: Angle.h:71
int const x0
Definition: saturated.cc:45
double asDegrees() const
Definition: Angle.h:124
virtual void shiftReferencePixel(double dx, double dy)
Move the pixel reference position by (dx, dy)
Definition: Wcs.cc:1161
boost::shared_ptr< RecordT > addNew()
Create a new record, add it to the end of the catalog, and return a pointer to it.
Definition: Catalog.h:469
virtual geom::Point2D skyToPixelImpl(geom::Angle sky1, geom::Angle sky2) const
Definition: Wcs.cc:738
table::Key< table::Point< double > > crval
Definition: Wcs.cc:1041
table::Key< table::Point< int > > pixel
table::Key< table::Array< Kernel::Pixel > > image
Definition: FixedKernel.cc:117
int d
Definition: KDTree.cc:89
A base class for objects that can be persisted via afw::table::io Archive classes.
Definition: Persistable.h:74
Field< T >::Value get(Key< T > const &key) const
Return the value of a field for the given key.
Definition: BaseRecord.h:123
A 2D linear coordinate transformation.
lsst::daf::base::PropertySet PropertySet
Definition: Wcs.cc:57
virtual std::string getPythonModule() const
Return the fully-qualified Python module that should be imported to guarantee that its factory is reg...
Definition: Wcs.cc:1082
virtual boost::shared_ptr< lsst::daf::base::PropertyList > getFitsMetadata() const
Return a PropertyList containing FITS header keywords that can be used to save the Wcs...
Definition: Wcs.cc:674
virtual Point2D reverseTransform(Point2D const &pixel) const
Definition: Wcs.cc:1304
AngleUnit const degrees
Definition: Angle.h:92
struct wcsprm * _wcsInfo
Definition: Wcs.h:390
int contains(Schema const &other, int flags=EQUAL_KEYS) const
Test whether the given schema is a subset of this.
virtual bool isPersistable() const
Whether the Wcs is persistable using afw::table::io archives.
Definition: Wcs.cc:1102
tbl::Schema schema
Definition: CoaddPsf.cc:324
An affine coordinate transformation consisting of a linear transformation and an offset.
virtual void rotateImageBy90(int nQuarter, lsst::afw::geom::Extent2I dimensions) const
Rotate image by nQuarter times 90 degrees.
Definition: Wcs.cc:623
table::Key< double > equinox
Definition: Wcs.cc:1046
bool operator==(Wcs const &other) const
Definition: Wcs.cc:496
#define CHECK_NULLS(a, b)
Definition: Wcs.cc:519
Eigen::Matrix2d getCDMatrix() const
Returns the CD matrix.
Definition: Wcs.cc:583
geom::Angle pixelScale() const
Returns the pixel scale [Angle/pixel].
Definition: Wcs.cc:731
table::Key< table::Point< double > > crpix
Definition: Wcs.cc:1042
const int STRLEN
Definition: Wcs.cc:65
lsst::afw::geom::Point2D getPixelOrigin() const
Returns CRPIX (corrected to LSST convention).
Definition: Wcs.cc:575
const int fitsToLsstPixels
Definition: Wcs.cc:79
Coord::Ptr makeCoord(CoordSystem const system, lsst::afw::geom::Angle const ra, lsst::afw::geom::Angle const dec, double const epoch)
Factory function to create a Coord of arbitrary type with decimal RA,Dec.
Definition: Coord.cc:1210
virtual void flipImage(int flipLR, int flipTB, lsst::afw::geom::Extent2I dimensions) const
Flip CD matrix around the y-axis.
Definition: Wcs.cc:601
boost::shared_ptr< lsst::daf::base::PropertyList > createTrivialWcsAsPropertySet(std::string const &wcsName, int const x0=0, int const y0=0)
Definition: Wcs.cc:1184
lsst::afw::coord::Coord::Ptr getSkyOrigin() const
Returns CRVAL. This need not be the centre of the image.
Definition: Wcs.cc:570
virtual geom::AffineTransform linearizePixelToSkyInternal(geom::Point2D const &pix, coord::Coord const &coord, geom::AngleUnit skyUnit) const
Definition: Wcs.cc:957
void shift(Extent< T, N > const &offset)
Shift the point by the given offset.
Definition: Point.h:113
int x
lsst::afw::image::Wcs Wcs
Definition: Wcs.cc:59
void _initWcs()
Definition: Wcs.cc:117
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
A vector of catalogs used by Persistable.
Definition: CatalogVector.h:26
int _relax
Degree of permissiveness for wcspih (0 for strict); see wcshdr.h for details.
Definition: Wcs.h:392
bool isFlipped() const
Definition: Wcs.cc:678
Base class for all records.
Definition: BaseRecord.h:27
geom::LinearTransform getLinearTransform() const
Definition: Wcs.cc:1014
A class used as a handle to a particular field in a table.
Definition: fwd.h:44
Wcs()
Construct an invalid Wcs given no arguments.
Definition: Wcs.cc:87
int status
const double PixelZeroPos
FITS uses 1.0, SDSS uses 0.5, LSST uses 0.0 (http://dev.lsstcorp.org/trac/wiki/BottomLeftPixelProposa...
Definition: ImageUtils.h:42
Class for storing generic metadata.
Definition: PropertySet.h:82
Virtual base class for 2D transforms.
Definition: XYTransform.h:48
int getAsInt(std::string const &name) const
Definition: PropertySet.cc:333
XYTransformFromWcsPair(boost::shared_ptr< Wcs const > dst, boost::shared_ptr< Wcs const > src)
Definition: Wcs.cc:1281
double pixArea(lsst::afw::geom::Point2D pix00) const
Sky area covered by a pixel at position pix00 in units of square degrees.
Definition: Wcs.cc:698
virtual void remove(std::string const &name)
Definition: PropertySet.cc:754
geom::AffineTransform linearizeSkyToPixel(coord::Coord const &coord, geom::AngleUnit skyUnit=geom::degrees) const
Return the local linear approximation to Wcs::skyToPixel at a point given in sky coordinates.
Definition: Wcs.cc:987
afw::table::Key< double > b
lsst::afw::image::XYTransformFromWcsPair XYTransformFromWcsPair
Definition: Wcs.cc:62
table::Key< std::string > cunit1
Definition: Wcs.cc:1048
Point< double, 2 > Point2D
Definition: Point.h:286
#define LSST_ARCHIVE_ASSERT(EXPR)
An assertion macro used to validate the structure of an InputArchive.
Definition: Persistable.h:47
CoordSystem makeCoordEnum(std::string const system)
A utility function to get the enum value of a coordinate system from a string name.
Definition: Coord.cc:116
std::string formatFitsProperties(boost::shared_ptr< lsst::daf::base::PropertySet const > const &prop)
double getAsDouble(std::string const &name) const
Definition: PropertySet.cc:406
double const PI
The ratio of a circle&#39;s circumference to diameter.
Definition: Angle.h:18
virtual void write(OutputArchiveHandle &handle) const
Write the object to one or more catalogs.
Definition: Wcs.cc:1084
A multi-catalog archive object used to load table::io::Persistable objects.
Definition: InputArchive.h:28
Coord::Ptr convert(CoordSystem system) const
Convert to a specified Coord type.
Definition: Coord.cc:636
int stripWcsKeywords(boost::shared_ptr< lsst::daf::base::PropertySet > const &metadata, boost::shared_ptr< Wcs const > const &wcs)
Definition: Wcs.cc:1250
int const y0
Definition: saturated.cc:45
int _wcshdrCtrl
Controls messages to stderr from wcshdr (0 for none); see wcshdr.h for details.
Definition: Wcs.h:394
lsst::afw::coord::Coord::Ptr CoordPtr
Definition: Wcs.cc:61
ImageT val
Definition: CR.cc:154
virtual Point2D forwardTransform(Point2D const &pixel) const
virtuals for forward and reverse transforms
Definition: Wcs.cc:1292
void saveCatalog(BaseCatalog const &catalog)
Save a catalog in the archive.
boost::shared_ptr< Wcs const > _dst
Definition: Wcs.h:445
dictionary Point
Definition: __init__.py:38
Include files required for standard LSST Exception handling.
int countFitsHeaderCards(boost::shared_ptr< lsst::daf::base::PropertySet const > const &prop)
table::Key< std::string > ctype1
Definition: Wcs.cc:1044
lsst::afw::geom::Point2D getPosition(lsst::afw::geom::AngleUnit unit=lsst::afw::geom::degrees) const
Return our contents in a Point2D object.
Definition: Coord.cc:476
table::Key< std::string > cunit2
Definition: Wcs.cc:1049
Functions to handle coordinates.