lsst.jointcal  16.0-23-gcb65559
ConstrainedPhotometryModel.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 <map>
26 #include <limits>
27 #include <vector>
28 #include <string>
29 
30 #include "lsst/log/Log.h"
31 
32 #include "astshim.h"
33 #include "astshim/ChebyMap.h"
37 #include "lsst/jointcal/CcdImage.h"
40 
41 namespace lsst {
42 namespace jointcal {
43 
44 unsigned ConstrainedPhotometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
45  unsigned index = firstIndex;
46  if (whatToFit.find("Model") == std::string::npos) {
47  LOGLS_WARN(_log, "assignIndices was called and Model is *not* in whatToFit");
48  return index;
49  }
50 
51  // If we got here, "Model" is definitely in whatToFit.
52  _fittingChips = (whatToFit.find("ModelChip") != std::string::npos);
53  _fittingVisits = (whatToFit.find("ModelVisit") != std::string::npos);
54  // If nothing more than "Model" is specified, it means fit everything.
55  if ((!_fittingChips) && (!_fittingVisits)) {
56  _fittingChips = _fittingVisits = true;
57  }
58 
59  if (_fittingChips) {
60  for (auto &idMapping : _chipMap) {
61  auto mapping = idMapping.second.get();
62  // Don't assign indices for fixed parameters.
63  if (mapping->isFixed()) continue;
64  mapping->setIndex(index);
65  index += mapping->getNpar();
66  }
67  }
68  if (_fittingVisits) {
69  for (auto &idMapping : _visitMap) {
70  auto mapping = idMapping.second.get();
71  mapping->setIndex(index);
72  index += mapping->getNpar();
73  }
74  }
75  for (auto &idMapping : _chipVisitMap) {
76  idMapping.second->setWhatToFit(_fittingChips, _fittingVisits);
77  }
78  return index;
79 }
80 
81 void ConstrainedPhotometryModel::offsetParams(Eigen::VectorXd const &delta) {
82  if (_fittingChips) {
83  for (auto &idMapping : _chipMap) {
84  auto mapping = idMapping.second.get();
85  // Don't offset indices for fixed parameters.
86  if (mapping->isFixed()) continue;
87  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
88  }
89  }
90  if (_fittingVisits) {
91  for (auto &idMapping : _visitMap) {
92  auto mapping = idMapping.second.get();
93  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
94  }
95  }
96 }
97 
99  for (auto &idMapping : _chipMap) {
100  idMapping.second.get()->freezeErrorTransform();
101  }
102  for (auto &idMapping : _visitMap) {
103  idMapping.second.get()->freezeErrorTransform();
104  }
105 }
106 
108  std::vector<unsigned> &indices) const {
109  auto mapping = findMapping(ccdImage);
110  mapping->getMappingIndices(indices);
111 }
112 
114  int total = 0;
115  for (auto &idMapping : _chipMap) {
116  total += idMapping.second->getNpar();
117  }
118  for (auto &idMapping : _visitMap) {
119  total += idMapping.second->getNpar();
120  }
121  return total;
122 }
123 
125  CcdImage const &ccdImage,
126  Eigen::VectorXd &derivatives) const {
127  auto mapping = findMapping(ccdImage);
128  mapping->computeParameterDerivatives(measuredStar, measuredStar.getInstFlux(), derivatives);
129 }
130 
131 namespace {
132 // Convert photoTransfo's way of storing Chebyshev coefficients into the format wanted by ChebyMap.
133 ndarray::Array<double, 2, 2> toChebyMapCoeffs(std::shared_ptr<PhotometryTransfoChebyshev> transfo) {
134  auto coeffs = transfo->getCoefficients();
135  // 4 x nPar: ChebyMap wants rows that look like (a_ij, 1, i, j) for out += a_ij*T_i(x)*T_j(y)
136  ndarray::Array<double, 2, 2> chebyCoeffs = allocate(ndarray::makeVector(transfo->getNpar(), 4));
137  Eigen::VectorXd::Index k = 0;
138  auto order = transfo->getOrder();
139  for (ndarray::Size j = 0; j <= order; ++j) {
140  ndarray::Size const iMax = order - j; // to save re-computing `i+j <= order` every inner step.
141  for (ndarray::Size i = 0; i <= iMax; ++i, ++k) {
142  chebyCoeffs[k][0] = coeffs[j][i];
143  chebyCoeffs[k][1] = 1;
144  chebyCoeffs[k][2] = i;
145  chebyCoeffs[k][3] = j;
146  }
147  }
148  return chebyCoeffs;
149 }
150 } // namespace
151 
153  for (auto &idMapping : _chipMap) {
154  idMapping.second->dump(stream);
155  stream << std::endl;
156  }
157  stream << std::endl;
158  for (auto &idMapping : _visitMap) {
159  idMapping.second->dump(stream);
160  stream << std::endl;
161  }
162 }
163 
165  auto idMapping = _chipVisitMap.find(ccdImage.getHashKey());
166  if (idMapping == _chipVisitMap.end())
168  "ConstrainedPhotometryModel cannot find CcdImage " + ccdImage.getName());
169  return idMapping->second.get();
170 }
171 
172 template <class ChipTransfo, class VisitTransfo, class ChipVisitMapping>
174  afw::geom::Box2D const &focalPlaneBBox, int visitOrder) {
175  // keep track of which chip we want to constrain (the one closest to the middle of the focal plane)
176  double minRadius2 = std::numeric_limits<double>::infinity();
177  CcdIdType constrainedChip = -1;
178 
179  // First initialize all visit and ccd transfos, before we make the ccdImage mappings.
180  for (auto const &ccdImage : ccdImageList) {
181  auto visit = ccdImage->getVisit();
182  auto chip = ccdImage->getCcdId();
183  auto visitPair = _visitMap.find(visit);
184  auto chipPair = _chipMap.find(chip);
185 
186  // If the chip is not in the map, add it, otherwise continue.
187  if (chipPair == _chipMap.end()) {
188  auto center = ccdImage->getDetector()->getCenter(afw::cameraGeom::FOCAL_PLANE);
189  double radius2 = std::pow(center.getX(), 2) + std::pow(center.getY(), 2);
190  if (radius2 < minRadius2) {
191  minRadius2 = radius2;
192  constrainedChip = chip;
193  }
194  auto photoCalib = ccdImage->getPhotoCalib();
195  // Use the single-frame processing calibration from the PhotoCalib as the default.
196  auto chipTransfo = std::make_unique<ChipTransfo>(initialChipCalibration(photoCalib));
197  _chipMap[chip] = std::make_shared<PhotometryMapping>(std::move(chipTransfo));
198  }
199  // If the visit is not in the map, add it, otherwise continue.
200  if (visitPair == _visitMap.end()) {
201  auto visitTransfo = std::make_unique<VisitTransfo>(visitOrder, focalPlaneBBox);
202  _visitMap[visit] = std::make_shared<PhotometryMapping>(std::move(visitTransfo));
203  }
204  }
205 
206  // Fix one chip mapping, to remove the degeneracy from the system.
207  _chipMap.at(constrainedChip)->setFixed(true);
208 
209  // Now create the ccdImage mappings, which are combinations of the chip/visit mappings above.
210  for (auto const &ccdImage : ccdImageList) {
211  auto visit = ccdImage->getVisit();
212  auto chip = ccdImage->getCcdId();
213  _chipVisitMap.emplace(ccdImage->getHashKey(),
214  std::make_unique<ChipVisitMapping>(_chipMap[chip], _visitMap[visit]));
215  }
216  LOGLS_INFO(_log, "Got " << _chipMap.size() << " chip mappings and " << _visitMap.size()
217  << " visit mappings; holding chip " << constrainedChip << " fixed ("
218  << getTotalParameters() << " total parameters).");
219  LOGLS_DEBUG(_log, "CcdImage map has " << _chipVisitMap.size() << " mappings, with "
220  << _chipVisitMap.bucket_count() << " buckets and a load factor of "
222 }
223 
224 // ConstrainedFluxModel methods
225 
227  MeasuredStar const &measuredStar) const {
228  return transform(ccdImage, measuredStar) - measuredStar.getFittedStar()->getFlux();
229 }
230 
231 double ConstrainedFluxModel::transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const {
232  auto mapping = findMapping(ccdImage);
233  return mapping->transform(measuredStar, measuredStar.getInstFlux());
234 }
235 
237  MeasuredStar const &measuredStar) const {
238  auto mapping = findMapping(ccdImage);
239  double tempErr = tweakFluxError(measuredStar);
240  return mapping->transformError(measuredStar, measuredStar.getInstFlux(), tempErr);
241 }
242 
244  auto oldPhotoCalib = ccdImage.getPhotoCalib();
245  auto detector = ccdImage.getDetector();
246  auto ccdBBox = detector->getBBox();
247  ChipVisitPhotometryMapping *mapping = dynamic_cast<ChipVisitPhotometryMapping *>(findMapping(ccdImage));
248  // There should be no way in which we can get to this point and not have a ChipVisitMapping,
249  // so blow up if we don't.
250  assert(mapping != nullptr);
251  auto pixToFocal = detector->getTransform(afw::cameraGeom::PIXELS, afw::cameraGeom::FOCAL_PLANE);
252  // We know it's a Chebyshev transfo because we created it as such, so blow up if it's not.
253  auto visitTransfo =
255  assert(visitTransfo != nullptr);
256  auto focalBBox = visitTransfo->getBBox();
257 
258  // Unravel our chebyshev coefficients to build an astshim::ChebyMap.
259  auto coeff_f = toChebyMapCoeffs(
260  std::dynamic_pointer_cast<PhotometryTransfoChebyshev>(mapping->getVisitMapping()->getTransfo()));
261  // Bounds are the bbox
262  std::vector<double> lowerBound = {focalBBox.getMinX(), focalBBox.getMinY()};
263  std::vector<double> upperBound = {focalBBox.getMaxX(), focalBBox.getMaxY()};
264 
265  afw::geom::TransformPoint2ToGeneric chebyTransform(ast::ChebyMap(coeff_f, 1, lowerBound, upperBound));
266 
267  // The chip part is easy: zoom map with the single value as the "zoom" factor.
269  ast::ZoomMap(1, mapping->getChipMapping()->getParameters()[0]));
270 
271  // Now stitch them all together.
272  auto transform = pixToFocal->then(chebyTransform)->then(zoomTransform);
273  // NOTE: TransformBoundedField does not yet implement mean(), so we have to compute it here.
274  // TODO: restore this calculation as part of DM-16305
275  // double mean = mapping->getChipMapping()->getParameters()[0] * visitTransfo->mean(ccdBBoxInFocal);
276  auto boundedField = std::make_shared<afw::math::TransformBoundedField>(ccdBBox, *transform);
277  return std::make_shared<afw::image::PhotoCalib>(oldPhotoCalib->getCalibrationMean(),
278  oldPhotoCalib->getCalibrationErr(), boundedField, false);
279 }
280 
281 // ConstrainedMagnitudeModel methods
282 
284  MeasuredStar const &measuredStar) const {
285  return transform(ccdImage, measuredStar) - measuredStar.getFittedStar()->getMag();
286 }
287 
289  MeasuredStar const &measuredStar) const {
290  auto mapping = findMapping(ccdImage);
291  return mapping->transform(measuredStar, measuredStar.getInstMag());
292 }
293 
295  MeasuredStar const &measuredStar) const {
296  auto mapping = findMapping(ccdImage);
297  double tempErr = tweakFluxError(measuredStar);
298  return mapping->transformError(measuredStar, measuredStar.getInstFlux(), tempErr);
299 }
300 
302  CcdImage const &ccdImage) const {
303  auto oldPhotoCalib = ccdImage.getPhotoCalib();
304  auto detector = ccdImage.getDetector();
305  auto ccdBBox = detector->getBBox();
306  ChipVisitPhotometryMapping *mapping = dynamic_cast<ChipVisitPhotometryMapping *>(findMapping(ccdImage));
307  // There should be no way in which we can get to this point and not have a ChipVisitMapping,
308  // so blow up if we don't.
309  assert(mapping != nullptr);
310  auto pixToFocal = detector->getTransform(afw::cameraGeom::PIXELS, afw::cameraGeom::FOCAL_PLANE);
311  // We know it's a Chebyshev transfo because we created it as such, so blow up if it's not.
312  auto visitTransfo =
314  assert(visitTransfo != nullptr);
315  auto focalBBox = visitTransfo->getBBox();
316 
317  // Unravel our chebyshev coefficients to build an astshim::ChebyMap.
318  auto coeff_f = toChebyMapCoeffs(
319  std::dynamic_pointer_cast<PhotometryTransfoChebyshev>(mapping->getVisitMapping()->getTransfo()));
320  // Bounds are the bbox
321  std::vector<double> lowerBound = {focalBBox.getMinX(), focalBBox.getMinY()};
322  std::vector<double> upperBound = {focalBBox.getMaxX(), focalBBox.getMaxY()};
323 
324  afw::geom::TransformPoint2ToGeneric chebyTransform(ast::ChebyMap(coeff_f, 1, lowerBound, upperBound));
325 
326  using namespace std::string_literals; // for operator""s to convert string literal->std::string
328  ast::MathMap(1, 1, {"y=pow(10.0,x/-2.5)"s}, {"x=-2.5*log10(y)"s}));
329 
330  // The chip part is easy: zoom map with the value (converted to a flux) as the "zoom" factor.
331  double calibrationMean = std::pow(10, mapping->getChipMapping()->getParameters()[0] / -2.5);
333  ast::ZoomMap(1, calibrationMean));
334 
335  // Now stitch them all together.
336  auto transform = pixToFocal->then(chebyTransform)->then(logTransform)->then(zoomTransform);
337  // NOTE: TransformBoundedField does not yet implement mean(), so we have to compute it here.
338  // TODO: restore this calculation as part of DM-16305
339  // double mean = calibrationMean * visitTransfo->mean(ccdBBoxInFocal);
340  auto boundedField = std::make_shared<afw::math::TransformBoundedField>(ccdBBox, *transform);
341  return std::make_shared<afw::image::PhotoCalib>(oldPhotoCalib->getCalibrationMean(),
342  oldPhotoCalib->getCalibrationErr(), boundedField, false);
343 }
344 
345 // explicit instantiation of templated function, so pybind11 can
348  afw::geom::Box2D const &, int);
351  CcdImageList const &, afw::geom::Box2D const &, int);
352 
353 } // namespace jointcal
354 } // namespace lsst
Photometric offset independent of position, defined as -2.5 * log(flux / fluxMag0).
Relates transfo(s) to their position in the fitting matrix and allows interaction with the transfo(s)...
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:79
CameraSysPrefix const PIXELS
std::shared_ptr< PhotometryMapping > getChipMapping() const
table::Key< int > detector
T endl(T... args)
T bucket_count(T... args)
T end(T... args)
nth-order 2d Chebyshev photometry transfo, times the input flux.
T load_factor(T... args)
Photometric offset independent of position, defined as (fluxMag0)^-1.
table::Key< double > calibrationMean
#define LOGLS_INFO(logger, message)
STL class.
unsigned assignIndices(std::string const &whatToFit, unsigned firstIndex) override
Assign indices in the full matrix to the parameters being fit in the mappings, starting at firstIndex...
T at(T... args)
std::shared_ptr< PhotometryMapping > getVisitMapping() const
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
int getTotalParameters() const override
Return the total number of parameters in this model.
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
std::shared_ptr< afw::image::PhotoCalib > getPhotoCalib() const
Return the exposure&#39;s photometric calibration.
Definition: CcdImage.h:160
Class for a simple mapping implementing a generic Gtransfo.
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return a pointer to the mapping associated with this ccdImage.
objects measured on actual images.
Definition: MeasuredStar.h:42
T dynamic_pointer_cast(T... args)
void initialize(CcdImageList const &ccdImageList, afw::geom::Box2D const &focalPlaneBBox, int visitOrder)
Initialize the chip, visit, and chipVisit mappings by creating appropriate transfos and mappings...
T infinity(T... args)
void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const override
Compute the parametric derivatives of this model.
T move(T... args)
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
T find(T... args)
T size(T... args)
void dump(std::ostream &stream=std::cout) const override
Dump the contents of the transfos, for debugging.
#define LSST_EXCEPT(type,...)
double tweakFluxError(jointcal::MeasuredStar const &measuredStar) const
Add a fraction of the instrumental flux to the instrumental flux error, in quadrature.
LOG_LOGGER _log
lsst.logging instance, to be created by a subclass so that messages have consistent name...
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the provided amounts (by -delta).
T pow(T... args)
#define LOGLS_DEBUG(logger, message)
T emplace(T... args)
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
nth-order 2d Chebyshev photometry transfo.
CcdImageKey getHashKey() const
Definition: CcdImage.h:151
virtual double initialChipCalibration(std::shared_ptr< afw::image::PhotoCalib const > photoCalib)=0
Return the initial calibration to use from this photoCalib.
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
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
nth-order 2d Chebyshev photometry transfo, plus the input flux.
std::shared_ptr< FittedStar > getFittedStar() const
Definition: MeasuredStar.h:108
CameraSys const FOCAL_PLANE
STL class.
void freezeErrorTransform() override
Once this routine has been called, the error transform is not modified by offsetParams().
virtual double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const =0
Return the on-sky transformed flux for measuredStar on ccdImage.
void getMappingIndices(CcdImage const &ccdImage, std::vector< unsigned > &indices) const override
Get how this set of parameters (of length Npar()) map into the "grand" fit.
A two-level photometric transform: one for the ccd and one for the visit.
#define LOGLS_WARN(logger, message)
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.