lsst.jointcal  14.0-17-ge3cc87b+9
CcdImage.cc
Go to the documentation of this file.
1 #include <assert.h>
2 #include <string>
3 #include <sstream>
4 #include <math.h>
5 
7 #include "lsst/pex/exceptions.h"
8 #include "lsst/afw/image/Image.h"
10 #include "lsst/afw/geom/Angle.h"
11 #include "lsst/afw/geom/Point.h"
12 
13 #include "lsst/log/Log.h"
14 #include "lsst/jointcal/CcdImage.h"
17 #include "lsst/jointcal/Gtransfo.h"
18 #include "lsst/jointcal/Point.h"
19 
20 namespace jointcal = lsst::jointcal;
21 namespace afwImg = lsst::afw::image;
22 
23 namespace {
24 LOG_LOGGER _log = LOG_GET("jointcal.CcdImage");
25 }
26 
27 namespace lsst {
28 namespace jointcal {
29 
30 void CcdImage::loadCatalog(afw::table::SourceCatalog const &catalog, std::string const &fluxField) {
31  auto xKey = catalog.getSchema().find<double>("slot_Centroid_x").key;
32  auto yKey = catalog.getSchema().find<double>("slot_Centroid_y").key;
33  auto xsKey = catalog.getSchema().find<float>("slot_Centroid_xSigma").key;
34  auto ysKey = catalog.getSchema().find<float>("slot_Centroid_ySigma").key;
35  auto mxxKey = catalog.getSchema().find<double>("slot_Shape_xx").key;
36  auto myyKey = catalog.getSchema().find<double>("slot_Shape_yy").key;
37  auto mxyKey = catalog.getSchema().find<double>("slot_Shape_xy").key;
38  auto fluxKey = catalog.getSchema().find<double>(fluxField + "_flux").key;
39  auto fluxErrKey = catalog.getSchema().find<double>(fluxField + "_fluxSigma").key;
40 
41  auto transform = _detector->getTransform(afw::cameraGeom::PIXELS, afw::cameraGeom::FOCAL_PLANE);
42 
43  _wholeCatalog.clear();
44  for (auto const &record : catalog) {
45  auto ms = std::make_shared<MeasuredStar>();
46  ms->x = record.get(xKey);
47  ms->y = record.get(yKey);
48  ms->vx = std::pow(record.get(xsKey), 2);
49  ms->vy = std::pow(record.get(ysKey), 2);
50  auto pointFocal = transform->applyForward(record.getCentroid());
51  ms->setXFocal(pointFocal.getX());
52  ms->setYFocal(pointFocal.getY());
53  /* the xy covariance is not provided in the input catalog: we
54  cook it up from the x and y position variance and the shape
55  measurements: */
56  double mxx = record.get(mxxKey);
57  double myy = record.get(myyKey);
58  double mxy = record.get(mxyKey);
59  ms->vxy = mxy * (ms->vx + ms->vy) / (mxx + myy);
60  if (ms->vx < 0 || ms->vy < 0 || (ms->vxy * ms->vxy) > (ms->vx * ms->vy)) {
61  LOGLS_WARN(_log, "Bad source detected in loadCatalog : " << ms->vx << " " << ms->vy << " "
62  << ms->vxy * ms->vxy << " "
63  << ms->vx * ms->vy);
64  continue;
65  }
66  ms->setId(record.getId());
67  ms->setInstFlux(record.get(fluxKey));
68  ms->setInstFluxErr(record.get(fluxErrKey));
69  // TODO: the below lines will be less clumsy once DM-4044 is cleaned up and we can say:
70  // TODO: instFluxToMaggies(ms->getInstFlux(), ms) (because ms will be derived from afw::geom::Point).
71  afw::geom::Point<double, 2> point(ms->x, ms->y);
72  auto flux = _photoCalib->instFluxToMaggies(ms->getInstFlux(), ms->getInstFluxErr(), point);
73  ms->setFlux(flux.value);
74  ms->setFluxErr(flux.err);
75  ms->mag = _photoCalib->instFluxToMagnitude(ms->getInstFlux(), point);
76  ms->setCcdImage(this);
77  _wholeCatalog.push_back(std::move(ms));
78  }
79  _wholeCatalog.setCcdImage(this);
80 }
81 
84  std::string const &filter, std::shared_ptr<afw::image::PhotoCalib> photoCalib,
85  std::shared_ptr<afw::cameraGeom::Detector> detector, int visit, int ccdId,
86  std::string const &fluxField)
87  : _ccdId(ccdId), _visit(visit), _photoCalib(photoCalib), _detector(detector), _filter(filter) {
88  loadCatalog(catalog, fluxField);
89 
90  Point lowerLeft(bbox.getMinX(), bbox.getMinY());
91  Point upperRight(bbox.getMaxX(), bbox.getMaxY());
92  _imageFrame = Frame(lowerLeft, upperRight);
93 
94  _readWcs.reset(new jointcal::TanSipPix2RaDec(jointcal::convertTanWcs(wcs)));
95  _inverseReadWcs = _readWcs->inverseTransfo(0.01, _imageFrame);
96 
98  out << visit << "_" << ccdId;
99  _name = out.str();
100 
101  _boresightRaDec = visitInfo->getBoresightRaDec();
102  _airMass = visitInfo->getBoresightAirmass();
103  _mjd = visitInfo->getDate().get(lsst::daf::base::DateTime::MJD);
104  double latitude = visitInfo->getObservatory().getLatitude();
105  _lstObs = visitInfo->getEra();
106  _hourAngle = visitInfo->getBoresightHourAngle();
107 
108  // lsstSim doesn't manage ERA (and thus Hour Angle) properly, so it's going to be NaN.
109  // Because we need the refraction vector later, go with 0 HA to prevent crashes on that NaN.
110  if (std::isnan(_hourAngle) == true) {
111  _hourAngle = 0;
112  }
113 
114  if (_airMass == 1)
115  _sineta = _coseta = _tgz = 0;
116  else {
117  double cosz = 1. / _airMass;
118  double sinz = sqrt(1 - cosz * cosz); // astronomers usually observe above the horizon
119  _tgz = sinz / cosz;
120  // TODO: as part of DM-12473, we can remove all of this and just call _visitInfo.getParallacticAngle()
121  double dec = _boresightRaDec.getLatitude();
122  // x/y components of refraction angle, eta.]
123  double yEta = sin(_hourAngle);
124  double xEta = cos(dec) * tan(latitude) - sin(dec) * cos(_hourAngle);
125  double eta = atan2(yEta, xEta);
126  _sineta = sin(eta);
127  _coseta = cos(eta);
128  }
129 }
130 
131 void CcdImage::setCommonTangentPoint(Point const &commonTangentPoint) {
132  _commonTangentPoint = commonTangentPoint;
133 
134  // use some other variable in case we later have to actually convert the
135  // as-read wcs:
136  const BaseTanWcs *tanWcs = _readWcs.get();
137 
138  /* we don't assume here that we know the internals of TanPix2RaDec:
139  to construct pix->TP, we do pix->sky->TP, although pix->sky
140  actually goes through TP */
141  GtransfoLin identity;
142  TanRaDec2Pix raDec2TP(identity, tanWcs->getTangentPoint());
143  _pix2TP = gtransfoCompose(&raDec2TP, tanWcs);
144  TanPix2RaDec CTP2RaDec(identity, commonTangentPoint);
145  _CTP2TP = gtransfoCompose(&raDec2TP, &CTP2RaDec);
146 
147  // jump from one TP to an other:
148  TanRaDec2Pix raDec2CTP(identity, commonTangentPoint);
149  TanPix2RaDec TP2RaDec(identity, tanWcs->getTangentPoint());
150  _TP2CTP = gtransfoCompose(&raDec2CTP, &TP2RaDec);
151  _sky2TP.reset(new TanRaDec2Pix(identity, tanWcs->getTangentPoint()));
152 
153  // this one is needed for matches :
154  _pix2CommonTangentPlane = gtransfoCompose(&raDec2CTP, tanWcs);
155 }
156 } // namespace jointcal
157 } // namespace lsst
#define LOGLS_WARN(logger, message)
Implements the (forward) SIP distorsion scheme.
Definition: Gtransfo.h:483
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:292
std::unique_ptr< Gtransfo > gtransfoCompose(const Gtransfo *left, const Gtransfo *right)
Returns a pointer to a composition.
Definition: Gtransfo.cc:367
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:443
table::Key< geom::Angle > latitude
lsst::afw::geom::Angle getLatitude() const
T atan2(T... args)
STL class.
T sin(T... args)
rectangle with sides parallel to axes.
Definition: Frame.h:19
Class for a simple mapping implementing a generic Gtransfo.
T str(T... args)
T cos(T... args)
T move(T... args)
T get(T... args)
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:517
T pow(T... args)
T tan(T... args)
Point getTangentPoint() const
The tangent point (in degrees)
Definition: Gtransfo.cc:1265
T dec(T... args)
T isnan(T... args)
T transform(T... args)
T sqrt(T... args)
#define LOG_GET(logger)
TanSipPix2RaDec convertTanWcs(const std::shared_ptr< lsst::afw::image::TanWcs > wcs)
Transform an afw TanWcs into a Gtransfo.
void setCommonTangentPoint(Point const &commonTangentPoint)
Sets the common tangent point and computes necessary transforms.
Definition: CcdImage.cc:131