lsst.jointcal  master-gc935ebf72c
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
CcdImage.cc
Go to the documentation of this file.
1 #include <assert.h>
2 #include <string>
3 #include <sstream>
4 #include <math.h>
5 
6 #include "lsst/log/Log.h"
10 #include "lsst/afw/image/Image.h"
11 #include "lsst/jointcal/Gtransfo.h"
12 #include "lsst/jointcal/Point.h"
13 #include "lsst/pex/exceptions.h"
14 #include "lsst/afw/geom/Angle.h"
15 
16 namespace jointcal = lsst::jointcal;
17 namespace afwImg = lsst::afw::image;
18 
19 namespace {
20 LOG_LOGGER _log = LOG_GET("jointcal.CcdImage");
21 }
22 
23 namespace lsst {
24 namespace jointcal {
25 
26 static double sq(double x) { return x * x; }
27 
28 void CcdImage::LoadCatalog(const lsst::afw::table::SortedCatalogT<lsst::afw::table::SourceRecord> &catalog,
29  const std::string &fluxField) {
30  auto xKey = catalog.getSchema().find<double>("slot_Centroid_x").key;
31  auto yKey = catalog.getSchema().find<double>("slot_Centroid_y").key;
32  auto xsKey = catalog.getSchema().find<float>("slot_Centroid_xSigma").key;
33  auto ysKey = catalog.getSchema().find<float>("slot_Centroid_ySigma").key;
34  auto mxxKey = catalog.getSchema().find<double>("slot_Shape_xx").key;
35  auto myyKey = catalog.getSchema().find<double>("slot_Shape_yy").key;
36  auto mxyKey = catalog.getSchema().find<double>("slot_Shape_xy").key;
37  auto fluxKey = catalog.getSchema().find<double>(fluxField + "_flux").key;
38  auto efluxKey = catalog.getSchema().find<double>(fluxField + "_fluxSigma").key;
39 
40  _wholeCatalog.clear();
41  for (auto const &record : catalog) {
42  auto ms = std::make_shared<MeasuredStar>();
43  ms->x = record.get(xKey);
44  ms->y = record.get(yKey);
45  ms->vx = sq(record.get(xsKey));
46  ms->vy = sq(record.get(ysKey));
47  /* the xy covariance is not provided in the input catalog: we
48  cook it up from the x and y position variance and the shape
49  measurements: */
50  double mxx = record.get(mxxKey);
51  double myy = record.get(myyKey);
52  double mxy = record.get(mxyKey);
53  ms->vxy = mxy * (ms->vx + ms->vy) / (mxx + myy);
54  if (ms->vx < 0 || ms->vy < 0 || (ms->vxy * ms->vxy) > (ms->vx * ms->vy)) {
55  LOGLS_WARN(_log, "Bad source detected in LoadCatalog : " << ms->vx << " " << ms->vy << " "
56  << ms->vxy * ms->vxy << " "
57  << ms->vx * ms->vy);
58  continue;
59  }
60  ms->setFlux(record.get(fluxKey));
61  ms->eflux = record.get(efluxKey);
62  ms->mag = _calib->getMagnitude(ms->getFlux());
63  ms->setCcdImage(this);
64  _wholeCatalog.push_back(std::move(ms));
65  }
66  _wholeCatalog.setCcdImage(this);
67 }
68 
69 CcdImage::CcdImage(lsst::afw::table::SortedCatalogT<lsst::afw::table::SourceRecord> &record,
70  const PTR(lsst::afw::image::TanWcs) wcs, const PTR(lsst::afw::image::VisitInfo) visitInfo,
71  const lsst::afw::geom::Box2I &bbox, const std::string &filter,
72  const PTR(lsst::afw::image::Calib) calib, const int &visit, const int &ccdId,
73  const std::string &fluxField)
74  : _ccdId(ccdId),
75  _visit(visit),
76  _calib(calib),
77  _filter(filter) {
78  LoadCatalog(record, fluxField);
79 
80  Point lowerLeft(bbox.getMinX(), bbox.getMinY());
81  Point upperRight(bbox.getMaxX(), bbox.getMaxY());
82  _imageFrame = Frame(lowerLeft, upperRight);
83 
84  _readWcs.reset(new jointcal::TanSipPix2RaDec(jointcal::convertTanWcs(wcs)));
85  _inverseReadWcs = _readWcs->inverseTransfo(0.01, _imageFrame);
86 
87  std::stringstream out;
88  out << visit << "_" << ccdId;
89  _name = out.str();
90 
91  _boresightRaDec = visitInfo->getBoresightRaDec();
92  _airMass = visitInfo->getBoresightAirmass();
93  _mjd = visitInfo->getDate().get(lsst::daf::base::DateTime::MJD);
94  double latitude = visitInfo->getObservatory().getLatitude();
95  _lstObs = visitInfo->getEra();
96  _hourAngle = visitInfo->getBoresightHourAngle();
97 
98  // lsstSim doesn't manage ERA (and thus Hour Angle) properly, so it's going to be NaN.
99  // Because we need the refraction vector later, go with 0 HA to prevent crashes on that NaN.
100  if (std::isnan(_hourAngle) == true) {
101  _hourAngle = 0;
102  }
103 
104  if (_airMass == 1)
105  _sineta = _coseta = _tgz = 0;
106  else {
107  double cosz = 1. / _airMass;
108  double sinz = sqrt(1 - cosz * cosz); // astronomers usually observe above the horizon
109  _tgz = sinz / cosz;
110  _sineta = cos(latitude) * sin(_hourAngle) / sinz;
111  _coseta = sqrt(1 - _sineta * _sineta);
112  if (_boresightRaDec.getDec() > latitude) _coseta = -_coseta;
113  }
114 }
115 
116 void CcdImage::setCommonTangentPoint(const Point &commonTangentPoint) {
117  _commonTangentPoint = commonTangentPoint;
118 
119  // use some other variable in case we later have to actually convert the
120  // as-read wcs:
121  const BaseTanWcs *tanWcs = _readWcs.get();
122 
123  /* we don't assume here that we know the internals of TanPix2RaDec:
124  to construct pix->TP, we do pix->sky->TP, although pix->sky
125  actually goes through TP */
126  GtransfoLin identity;
127  TanRaDec2Pix raDec2TP(identity, tanWcs->getTangentPoint());
128  _pix2TP = gtransfoCompose(&raDec2TP, tanWcs);
129  TanPix2RaDec CTP2RaDec(identity, commonTangentPoint);
130  _CTP2TP = gtransfoCompose(&raDec2TP, &CTP2RaDec);
131 
132  // jump from one TP to an other:
133  TanRaDec2Pix raDec2CTP(identity, commonTangentPoint);
134  TanPix2RaDec TP2RaDec(identity, tanWcs->getTangentPoint());
135  _TP2CTP = gtransfoCompose(&raDec2CTP, &TP2RaDec);
136  _sky2TP.reset(new TanRaDec2Pix(identity, tanWcs->getTangentPoint()));
137 
138  // this one is needed for matches :
139  _pix2CommonTangentPlane = gtransfoCompose(&raDec2CTP, tanWcs);
140 }
141 } // namespace jointcal
142 } // namespace lsst
Implements the (forward) SIP distorsion scheme.
Definition: Gtransfo.h:479
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:288
std::unique_ptr< Gtransfo > gtransfoCompose(const Gtransfo *left, const Gtransfo *right)
Returns a pointer to a composition.
Definition: Gtransfo.cc:367
void setCommonTangentPoint(const Point &commonTangentPoint)
Sets the common tangent point and computes necessary transforms.
Definition: CcdImage.cc:116
A point in a plane.
Definition: Point.h:13
the transformation that handles pix to sideral transfos (Gnomonic, possibly with polynomial distortio...
Definition: Gtransfo.h:439
CcdImage(lsst::afw::table::SortedCatalogT< lsst::afw::table::SourceRecord > &record, const boost::shared_ptr< lsst::afw::image::TanWcs > wcs, const boost::shared_ptr< lsst::afw::image::VisitInfo > visitInfo, const lsst::afw::geom::Box2I &bbox, const std::string &filter, const boost::shared_ptr< lsst::afw::image::Calib > calib, const int &visit, const int &ccd, const std::string &fluxField)
Definition: CcdImage.cc:69
rectangle with sides parallel to axes.
Definition: Frame.h:19
void setCcdImage(const CcdImage *_ccdImage)
Definition: MeasuredStar.cc:45
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:513
Point getTangentPoint() const
The tangent point (in degrees)
Definition: Gtransfo.cc:1266
TanSipPix2RaDec convertTanWcs(const std::shared_ptr< lsst::afw::image::TanWcs > wcs)
Transform an afw TanWcs into a Gtransfo.