6 #include "lsst/log/Log.h"
10 #include "lsst/afw/image/Image.h"
13 #include "lsst/pex/exceptions.h"
14 #include "lsst/afw/geom/Angle.h"
16 namespace jointcal = lsst::jointcal;
17 namespace afwImg = lsst::afw::image;
20 LOG_LOGGER _log = LOG_GET(
"jointcal.CcdImage");
26 static double sq(
double x) {
return x * x; }
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;
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));
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 <<
" "
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));
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)
78 LoadCatalog(record, fluxField);
80 Point lowerLeft(bbox.getMinX(), bbox.getMinY());
81 Point upperRight(bbox.getMaxX(), bbox.getMaxY());
82 _imageFrame =
Frame(lowerLeft, upperRight);
85 _inverseReadWcs = _readWcs->inverseTransfo(0.01, _imageFrame);
87 std::stringstream out;
88 out << visit <<
"_" << ccdId;
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();
100 if (std::isnan(_hourAngle) ==
true) {
105 _sineta = _coseta = _tgz = 0;
107 double cosz = 1. / _airMass;
108 double sinz = sqrt(1 - cosz * cosz);
110 _sineta = cos(latitude) * sin(_hourAngle) / sinz;
111 _coseta = sqrt(1 - _sineta * _sineta);
112 if (_boresightRaDec.getDec() > latitude) _coseta = -_coseta;
117 _commonTangentPoint = commonTangentPoint;
Implements the (forward) SIP distorsion scheme.
implements the linear transformations (6 real coefficients).
std::unique_ptr< Gtransfo > gtransfoCompose(const Gtransfo *left, const Gtransfo *right)
Returns a pointer to a composition.
void setCommonTangentPoint(const Point &commonTangentPoint)
Sets the common tangent point and computes necessary transforms.
the transformation that handles pix to sideral transfos (Gnomonic, possibly with polynomial distortio...
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)
rectangle with sides parallel to axes.
void setCcdImage(const CcdImage *_ccdImage)
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Point getTangentPoint() const
The tangent point (in degrees)
TanSipPix2RaDec convertTanWcs(const std::shared_ptr< lsst::afw::image::TanWcs > wcs)
Transform an afw TanWcs into a Gtransfo.