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