lsst.jointcal  16.0-15-g8e16a51+12
ConstrainedAstrometryModel.cc
Go to the documentation of this file.
1 #include "astshim.h"
2 #include "lsst/afw/geom.h"
3 #include "lsst/log/Log.h"
11 
12 #include "lsst/pex/exceptions.h"
14 
15 #include <memory>
16 #include <string>
17 #include <iostream>
18 
19 namespace {
20 LOG_LOGGER _log = LOG_GET("jointcal.ConstrainedAstrometryModel");
21 }
22 
23 namespace {
24 // Append the keys of this map into a comma-separated string.
25 template <typename KeyType, typename ValueType>
26 void outputMapKeys(std::map<KeyType, ValueType> const &map, std::ostream &os) {
27  bool first = true;
28  os << "[";
29  for (auto const &i : map) {
30  if (first)
31  first = false;
32  else
33  os << ", ";
34  os << i.first;
35  }
36  os << "]";
37 }
38 } // namespace
39 
40 namespace lsst {
41 namespace jointcal {
42 
43 ConstrainedAstrometryModel::ConstrainedAstrometryModel(
45  int chipOrder, int visitOrder)
46  : _sky2TP(projectionHandler) {
47  // keep track of which chip we want to hold fixed (the one closest to the middle of the focal plane)
48  double minRadius2 = std::numeric_limits<double>::infinity();
49  CcdIdType constrainedChip = -1;
50 
51  // first loop to initialize all visit and chip transfos.
52  for (auto &ccdImage : ccdImageList) {
53  const CcdImage &im = *ccdImage;
54  auto visit = im.getVisit();
55  auto chip = im.getCcdId();
56  auto visitp = _visitMap.find(visit);
57  if (visitp == _visitMap.end()) {
58  _visitMap[visit] = std::make_shared<SimplePolyMapping>(GtransfoLin(), GtransfoPoly(visitOrder));
59  }
60  auto chipp = _chipMap.find(chip);
61  if (chipp == _chipMap.end()) {
62  auto center = ccdImage->getDetector()->getCenter(afw::cameraGeom::FOCAL_PLANE);
63  double radius2 = std::pow(center.getX(), 2) + std::pow(center.getY(), 2);
64  if (radius2 < minRadius2) {
65  minRadius2 = radius2;
66  constrainedChip = chip;
67  }
68 
69  auto pixelsToFocal =
71  Frame const &frame = im.getImageFrame();
72  // construct the chip gtransfo by approximating the pixel->Focal afw::geom::Transform.
73  GtransfoPoly pol = GtransfoPoly(pixelsToFocal, frame, chipOrder);
74  GtransfoLin shiftAndNormalize = normalizeCoordinatesTransfo(frame);
75  _chipMap[chip] = std::make_shared<SimplePolyMapping>(shiftAndNormalize,
76  pol * shiftAndNormalize.inverted());
77  }
78  }
79 
80  // Hold the "central" chip map fixed and don't fit it, to remove a degeneracy.
81  _chipMap.at(constrainedChip)->setToBeFit(false);
82 
83  // now, second loop to set the mappings of the CCdImages
84  for (auto &ccdImage : ccdImageList) {
85  const CcdImage &im = *ccdImage;
86  auto visit = im.getVisit();
87  auto chip = im.getCcdId();
88 
89  // check that the chip_indexed part was indeed assigned
90  // (i.e. the reference visit was complete)
91  if (_chipMap.find(chip) == _chipMap.end()) {
92  LOGLS_WARN(_log, "Chip " << chip << " is missing in the reference exposure, expect troubles.");
94  _chipMap[chip] = std::make_shared<SimplePolyMapping>(norm, GtransfoPoly(chipOrder));
95  }
96  _mappings[ccdImage->getHashKey()] =
97  std::make_unique<TwoTransfoMapping>(_chipMap[chip], _visitMap[visit]);
98  }
99  LOGLS_INFO(_log, "Got " << _chipMap.size() << " chip mappings and " << _visitMap.size()
100  << " visit mappings; holding chip " << constrainedChip << " fixed ("
101  << getTotalParameters() << " total parameters).");
102  LOGLS_DEBUG(_log, "CcdImage map has " << _mappings.size() << " mappings, with "
103  << _mappings.bucket_count() << " buckets and a load factor of "
104  << _mappings.load_factor());
105 }
106 
108  return findMapping(ccdImage);
109 }
110 
115 unsigned ConstrainedAstrometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
116  unsigned index = firstIndex;
117  if (whatToFit.find("Distortions") == std::string::npos) {
118  LOGLS_ERROR(_log, "assignIndices was called and Distortions is *not* in whatToFit");
119  return 0;
120  }
121  // if we get here "Distortions" is in whatToFit
122  _fittingChips = (whatToFit.find("DistortionsChip") != std::string::npos);
123  _fittingVisits = (whatToFit.find("DistortionsVisit") != std::string::npos);
124  // If nothing more than "Distortions" is specified, it means all:
125  if ((!_fittingChips) && (!_fittingVisits)) {
126  _fittingChips = _fittingVisits = true;
127  }
128  if (_fittingChips)
129  for (auto &i : _chipMap) {
130  i.second->setIndex(index);
131  index += i.second->getNpar();
132  }
133  if (_fittingVisits)
134  for (auto &i : _visitMap) {
135  i.second->setIndex(index);
136  index += i.second->getNpar();
137  }
138  // Tell the mappings which derivatives they will have to fill:
139  for (auto &i : _mappings) {
140  i.second->setWhatToFit(_fittingChips, _fittingVisits);
141  }
142  return index;
143 }
144 
145 void ConstrainedAstrometryModel::offsetParams(Eigen::VectorXd const &delta) {
146  if (_fittingChips)
147  for (auto &i : _chipMap) {
148  auto mapping = i.second.get();
149  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
150  }
151  if (_fittingVisits)
152  for (auto &i : _visitMap) {
153  auto mapping = i.second.get();
154  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
155  }
156 }
157 
159  for (auto i = _visitMap.begin(); i != _visitMap.end(); ++i) i->second->freezeErrorTransform();
160  for (auto i = _chipMap.begin(); i != _chipMap.end(); ++i) i->second->freezeErrorTransform();
161 }
162 
164  auto chipp = _chipMap.find(chip);
165  if (chipp == _chipMap.end()) {
166  std::stringstream errMsg;
167  errMsg << "No such chipId: " << chip << " among ";
168  outputMapKeys(_chipMap, errMsg);
169  std::cout << std::endl;
170  throw pexExcept::InvalidParameterError(errMsg.str());
171  }
172  return chipp->second->getTransfo();
173 }
174 
175 // Array of visits involved in the solution.
178  res.reserve(_visitMap.size());
179  for (auto i = _visitMap.begin(); i != _visitMap.end(); ++i) res.push_back(i->first);
180  return res;
181 }
182 
184  auto visitp = _visitMap.find(visit);
185  if (visitp == _visitMap.end()) {
186  std::stringstream errMsg;
187  errMsg << "No such visitId: " << visit << " among ";
188  outputMapKeys(_visitMap, errMsg);
189  std::cout << std::endl;
190  throw pexExcept::InvalidParameterError(errMsg.str());
191  }
192  return visitp->second->getTransfo();
193 }
194 
196  int total = 0;
197  for (auto &i : _chipMap) {
198  total += i.second->getNpar();
199  }
200  for (auto &i : _visitMap) {
201  total += i.second->getNpar();
202  }
203  return total;
204 }
205 
207  auto proj = std::dynamic_pointer_cast<const TanRaDec2Pix>(getSky2TP(ccdImage));
208  jointcal::Point tangentPoint(proj->getTangentPoint());
209 
210  auto imageFrame = ccdImage.getImageFrame();
211  auto pixelsToFocal = getChipTransfo(ccdImage.getCcdId()).toAstMap(imageFrame);
212  jointcal::Frame focalBox = getChipTransfo(ccdImage.getCcdId()).apply(imageFrame, false);
213  auto focalToIwc = getVisitTransfo(ccdImage.getVisit()).toAstMap(focalBox);
214 
215  ast::Frame pixelFrame(2, "Domain=PIXELS");
216  ast::Frame focalFrame(2, "Domain=FOCAL");
217  ast::Frame iwcFrame(2, "Domain=IWC");
218 
219  // make a basic SkyWcs and extract the IWC portion
220  auto iwcToSkyWcs = afw::geom::makeSkyWcs(
221  afw::geom::Point2D(0, 0),
222  afw::geom::SpherePoint(tangentPoint.x, tangentPoint.y, afw::geom::degrees),
223  afw::geom::makeCdMatrix(1.0 * afw::geom::degrees, 0 * afw::geom::degrees, true));
224  auto iwcToSkyMap = iwcToSkyWcs->getFrameDict()->getMapping("PIXELS", "SKY");
225  auto skyFrame = iwcToSkyWcs->getFrameDict()->getFrame("SKY");
226 
227  ast::FrameDict frameDict(pixelFrame);
228  frameDict.addFrame("PIXELS", *pixelsToFocal, focalFrame);
229  frameDict.addFrame("FOCAL", *focalToIwc, iwcFrame);
230  frameDict.addFrame("IWC", *iwcToSkyMap, *skyFrame);
231  return std::make_shared<afw::geom::SkyWcs>(frameDict);
232 }
233 
234 AstrometryMapping *ConstrainedAstrometryModel::findMapping(CcdImage const &ccdImage) const {
235  auto i = _mappings.find(ccdImage.getHashKey());
236  if (i == _mappings.end())
238  "ConstrainedAstrometryModel cannot find CcdImage " + ccdImage.getName());
239  return i->second.get();
240 }
241 
242 } // namespace jointcal
243 } // namespace lsst
#define LOGLS_WARN(logger, message)
unsigned assignIndices(std::string const &whatToFit, unsigned firstIndex) override
Positions the various parameter sets into the parameter vector, starting at firstIndex.
VisitIdType getVisit() const
returns visit ID
Definition: CcdImage.h:124
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:391
A point in a plane.
Definition: Point.h:13
const std::shared_ptr< Gtransfo const > getSky2TP(CcdImage const &ccdImage) const override
The mapping of sky coordinates (i.e.
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:56
CameraSysPrefix const PIXELS
T endl(T... args)
T norm(const T &x)
std::vector< VisitIdType > getVisits() const
Access to array of visits involved in the solution.
CcdIdType getCcdId() const
returns ccd ID
Definition: CcdImage.h:121
Gtransfo const & getVisitTransfo(VisitIdType const &visit) const
Access to mappings.
Polynomial transformation class.
Definition: Gtransfo.h:253
GtransfoLin normalizeCoordinatesTransfo(const Frame &frame)
Returns the transformation that maps the input frame along both axes to [-1,1].
Definition: Gtransfo.cc:764
STL class.
AstrometryMapping const * getMapping(CcdImage const &) const override
Mapping associated to a given CcdImage.
STL class.
pairs of points
T push_back(T... args)
first
rectangle with sides parallel to axes.
Definition: Frame.h:15
#define LOGLS_DEBUG(logger, message)
Class for a simple mapping implementing a generic Gtransfo.
void addFrame(int iframe, Mapping const &map, Frame const &frame) override
std::shared_ptr< afw::geom::SkyWcs > makeSkyWcs(CcdImage const &ccdImage) const override
Make a SkyWcs that contains this model.
T str(T... args)
void freezeErrorTransform() override
From there on, measurement errors are propagated using the current transfos (and no longer evolve)...
T dynamic_pointer_cast(T... args)
T infinity(T... args)
Gtransfo const & getChipTransfo(CcdIdType const chip) const
Access to mappings.
#define LOGLS_INFO(logger, message)
T find(T... args)
#define LSST_EXCEPT(type,...)
STL class.
void offsetParams(Eigen::VectorXd const &Delta) override
Dispaches the offsets after a fit step into the actual locations of parameters.
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:653
int getTotalParameters() const override
Return the total number of parameters in this model.
T pow(T... args)
a virtual (interface) class for geometric transformations.
Definition: Gtransfo.h:42
std::shared_ptr< SkyWcs > makeSkyWcs(daf::base::PropertySet &metadata, bool strip=false)
Eigen::Matrix2d makeCdMatrix(lsst::geom::Angle const &scale, lsst::geom::Angle const &orientation=0 *lsst::geom::degrees, bool flipX=false)
int VisitIdType
Definition: CcdImage.h:25
virtual class needed in the abstraction of the distortion model
CcdImageKey getHashKey() const
Definition: CcdImage.h:128
std::shared_ptr< afw::cameraGeom::Detector > getDetector() const
Definition: CcdImage.h:126
Handler of an actual image from a single CCD.
Definition: CcdImage.h:41
#define LOG_GET(logger)
CameraSys const FOCAL_PLANE
STL class.
#define LOGLS_ERROR(logger, message)
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:166
T reserve(T... args)