lsst.jointcal  16.0-24-gd8faad8+2
ConstrainedAstrometryModel.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 /*
3  * This file is part of jointcal.
4  *
5  * Developed for the LSST Data Management System.
6  * This product includes software developed by the LSST Project
7  * (https://www.lsst.org).
8  * See the COPYRIGHT file at the top-level directory of this distribution
9  * for details of code ownership.
10  *
11  * This program is free software: you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License as published by
13  * the Free Software Foundation, either version 3 of the License, or
14  * (at your option) any later version.
15  *
16  * This program is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19  * GNU General Public License for more details.
20  *
21  * You should have received a copy of the GNU General Public License
22  * along with this program. If not, see <https://www.gnu.org/licenses/>.
23  */
24 
25 #include "astshim.h"
26 #include "lsst/afw/geom.h"
27 #include "lsst/log/Log.h"
30 #include "lsst/jointcal/CcdImage.h"
32 #include "lsst/jointcal/Gtransfo.h"
35 
36 #include "lsst/pex/exceptions.h"
38 
39 #include <memory>
40 #include <string>
41 #include <iostream>
42 
43 namespace {
44 LOG_LOGGER _log = LOG_GET("jointcal.ConstrainedAstrometryModel");
45 }
46 
47 namespace {
48 // Append the keys of this map into a comma-separated string.
49 template <typename KeyType, typename ValueType>
50 void outputMapKeys(std::map<KeyType, ValueType> const &map, std::ostream &os) {
51  bool first = true;
52  os << "[";
53  for (auto const &i : map) {
54  if (first)
55  first = false;
56  else
57  os << ", ";
58  os << i.first;
59  }
60  os << "]";
61 }
62 } // namespace
63 
64 namespace lsst {
65 namespace jointcal {
66 
67 ConstrainedAstrometryModel::ConstrainedAstrometryModel(
69  int chipOrder, int visitOrder)
70  : _sky2TP(projectionHandler) {
71  // keep track of which chip we want to hold fixed (the one closest to the middle of the focal plane)
72  double minRadius2 = std::numeric_limits<double>::infinity();
73  CcdIdType constrainedChip = -1;
74 
75  // first loop to initialize all visit and chip transfos.
76  for (auto &ccdImage : ccdImageList) {
77  const CcdImage &im = *ccdImage;
78  auto visit = im.getVisit();
79  auto chip = im.getCcdId();
80  auto visitp = _visitMap.find(visit);
81  if (visitp == _visitMap.end()) {
82  _visitMap[visit] = std::make_shared<SimplePolyMapping>(GtransfoLin(), GtransfoPoly(visitOrder));
83  }
84  auto chipp = _chipMap.find(chip);
85  if (chipp == _chipMap.end()) {
86  auto center = ccdImage->getDetector()->getCenter(afw::cameraGeom::FOCAL_PLANE);
87  double radius2 = std::pow(center.getX(), 2) + std::pow(center.getY(), 2);
88  if (radius2 < minRadius2) {
89  minRadius2 = radius2;
90  constrainedChip = chip;
91  }
92 
93  auto pixelsToFocal =
95  Frame const &frame = im.getImageFrame();
96  // construct the chip gtransfo by approximating the pixel->Focal afw::geom::Transform.
97  GtransfoPoly pol = GtransfoPoly(pixelsToFocal, frame, chipOrder);
98  GtransfoLin shiftAndNormalize = normalizeCoordinatesTransfo(frame);
99  _chipMap[chip] = std::make_shared<SimplePolyMapping>(shiftAndNormalize,
100  pol * shiftAndNormalize.inverted());
101  }
102  }
103 
104  // Hold the "central" chip map fixed and don't fit it, to remove a degeneracy.
105  _chipMap.at(constrainedChip)->setToBeFit(false);
106 
107  // now, second loop to set the mappings of the CCdImages
108  for (auto &ccdImage : ccdImageList) {
109  const CcdImage &im = *ccdImage;
110  auto visit = im.getVisit();
111  auto chip = im.getCcdId();
112 
113  // check that the chip_indexed part was indeed assigned
114  // (i.e. the reference visit was complete)
115  if (_chipMap.find(chip) == _chipMap.end()) {
116  LOGLS_WARN(_log, "Chip " << chip << " is missing in the reference exposure, expect troubles.");
118  _chipMap[chip] = std::make_shared<SimplePolyMapping>(norm, GtransfoPoly(chipOrder));
119  }
120  _mappings[ccdImage->getHashKey()] =
121  std::make_unique<TwoTransfoMapping>(_chipMap[chip], _visitMap[visit]);
122  }
123  LOGLS_INFO(_log, "Got " << _chipMap.size() << " chip mappings and " << _visitMap.size()
124  << " visit mappings; holding chip " << constrainedChip << " fixed ("
125  << getTotalParameters() << " total parameters).");
126  LOGLS_DEBUG(_log, "CcdImage map has " << _mappings.size() << " mappings, with "
127  << _mappings.bucket_count() << " buckets and a load factor of "
128  << _mappings.load_factor());
129 }
130 
132  return findMapping(ccdImage);
133 }
134 
139 unsigned ConstrainedAstrometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
140  unsigned index = firstIndex;
141  if (whatToFit.find("Distortions") == std::string::npos) {
142  LOGLS_ERROR(_log, "assignIndices was called and Distortions is *not* in whatToFit");
143  return 0;
144  }
145  // if we get here "Distortions" is in whatToFit
146  _fittingChips = (whatToFit.find("DistortionsChip") != std::string::npos);
147  _fittingVisits = (whatToFit.find("DistortionsVisit") != std::string::npos);
148  // If nothing more than "Distortions" is specified, it means all:
149  if ((!_fittingChips) && (!_fittingVisits)) {
150  _fittingChips = _fittingVisits = true;
151  }
152  if (_fittingChips)
153  for (auto &i : _chipMap) {
154  i.second->setIndex(index);
155  index += i.second->getNpar();
156  }
157  if (_fittingVisits)
158  for (auto &i : _visitMap) {
159  i.second->setIndex(index);
160  index += i.second->getNpar();
161  }
162  // Tell the mappings which derivatives they will have to fill:
163  for (auto &i : _mappings) {
164  i.second->setWhatToFit(_fittingChips, _fittingVisits);
165  }
166  return index;
167 }
168 
169 void ConstrainedAstrometryModel::offsetParams(Eigen::VectorXd const &delta) {
170  if (_fittingChips)
171  for (auto &i : _chipMap) {
172  auto mapping = i.second.get();
173  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
174  }
175  if (_fittingVisits)
176  for (auto &i : _visitMap) {
177  auto mapping = i.second.get();
178  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
179  }
180 }
181 
183  for (auto i = _visitMap.begin(); i != _visitMap.end(); ++i) i->second->freezeErrorTransform();
184  for (auto i = _chipMap.begin(); i != _chipMap.end(); ++i) i->second->freezeErrorTransform();
185 }
186 
188  auto chipp = _chipMap.find(chip);
189  if (chipp == _chipMap.end()) {
190  std::stringstream errMsg;
191  errMsg << "No such chipId: " << chip << " among ";
192  outputMapKeys(_chipMap, errMsg);
193  std::cout << std::endl;
194  throw pexExcept::InvalidParameterError(errMsg.str());
195  }
196  return chipp->second->getTransfo();
197 }
198 
199 // Array of visits involved in the solution.
202  res.reserve(_visitMap.size());
203  for (auto i = _visitMap.begin(); i != _visitMap.end(); ++i) res.push_back(i->first);
204  return res;
205 }
206 
208  auto visitp = _visitMap.find(visit);
209  if (visitp == _visitMap.end()) {
210  std::stringstream errMsg;
211  errMsg << "No such visitId: " << visit << " among ";
212  outputMapKeys(_visitMap, errMsg);
213  std::cout << std::endl;
214  throw pexExcept::InvalidParameterError(errMsg.str());
215  }
216  return visitp->second->getTransfo();
217 }
218 
220  int total = 0;
221  for (auto &i : _chipMap) {
222  total += i.second->getNpar();
223  }
224  for (auto &i : _visitMap) {
225  total += i.second->getNpar();
226  }
227  return total;
228 }
229 
231  auto proj = std::dynamic_pointer_cast<const TanRaDec2Pix>(getSky2TP(ccdImage));
232  jointcal::Point tangentPoint(proj->getTangentPoint());
233 
234  auto imageFrame = ccdImage.getImageFrame();
235  auto pixelsToFocal = getChipTransfo(ccdImage.getCcdId()).toAstMap(imageFrame);
236  jointcal::Frame focalBox = getChipTransfo(ccdImage.getCcdId()).apply(imageFrame, false);
237  auto focalToIwc = getVisitTransfo(ccdImage.getVisit()).toAstMap(focalBox);
238 
239  ast::Frame pixelFrame(2, "Domain=PIXELS");
240  ast::Frame focalFrame(2, "Domain=FOCAL");
241  ast::Frame iwcFrame(2, "Domain=IWC");
242 
243  // make a basic SkyWcs and extract the IWC portion
244  auto iwcToSkyWcs = afw::geom::makeSkyWcs(
245  afw::geom::Point2D(0, 0),
246  afw::geom::SpherePoint(tangentPoint.x, tangentPoint.y, afw::geom::degrees),
247  afw::geom::makeCdMatrix(1.0 * afw::geom::degrees, 0 * afw::geom::degrees, true));
248  auto iwcToSkyMap = iwcToSkyWcs->getFrameDict()->getMapping("PIXELS", "SKY");
249  auto skyFrame = iwcToSkyWcs->getFrameDict()->getFrame("SKY");
250 
251  ast::FrameDict frameDict(pixelFrame);
252  frameDict.addFrame("PIXELS", *pixelsToFocal, focalFrame);
253  frameDict.addFrame("FOCAL", *focalToIwc, iwcFrame);
254  frameDict.addFrame("IWC", *iwcToSkyMap, *skyFrame);
255  return std::make_shared<afw::geom::SkyWcs>(frameDict);
256 }
257 
258 AstrometryMapping *ConstrainedAstrometryModel::findMapping(CcdImage const &ccdImage) const {
259  auto i = _mappings.find(ccdImage.getHashKey());
260  if (i == _mappings.end())
262  "ConstrainedAstrometryModel cannot find CcdImage " + ccdImage.getName());
263  return i->second.get();
264 }
265 
266 } // namespace jointcal
267 } // namespace lsst
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:147
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:414
A point in a plane.
Definition: Point.h:36
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:79
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:144
Gtransfo const & getVisitTransfo(VisitIdType const &visit) const
Access to mappings.
Polynomial transformation class.
Definition: Gtransfo.h:276
GtransfoLin normalizeCoordinatesTransfo(const Frame &frame)
Returns the transformation that maps the input frame along both axes to [-1,1].
Definition: Gtransfo.cc:788
STL class.
AstrometryMapping const * getMapping(CcdImage const &) const override
Mapping associated to a given CcdImage.
#define LOGLS_INFO(logger, message)
STL class.
pairs of points
T push_back(T... args)
void addFrame(int iframe, Mapping const &map, Frame const &frame) override
first
rectangle with sides parallel to axes.
Definition: Frame.h:38
Class for a simple mapping implementing a generic Gtransfo.
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.
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:676
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:65
std::shared_ptr< SkyWcs > makeSkyWcs(daf::base::PropertySet &metadata, bool strip=false)
#define LOGLS_DEBUG(logger, message)
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:48
virtual class needed in the abstraction of the distortion model
CcdImageKey getHashKey() const
Definition: CcdImage.h:151
std::shared_ptr< afw::cameraGeom::Detector > getDetector() const
Definition: CcdImage.h:149
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
CameraSys const FOCAL_PLANE
STL class.
#define LOG_GET(logger)
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:189
#define LOGLS_ERROR(logger, message)
#define LOGLS_WARN(logger, message)
T reserve(T... args)