lsst.jointcal  16.0-27-ge118ca2
SimplePhotometryModel.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 <iostream>
26 #include <cmath>
27 
28 #include "lsst/log/Log.h"
32 #include "lsst/jointcal/CcdImage.h"
34 
35 namespace {
36 double instFluxFromMag(double mag) { return std::pow(10, mag / -2.5); }
37 } // namespace
38 
39 namespace lsst {
40 namespace jointcal {
41 
42 unsigned SimplePhotometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
43  unsigned ipar = firstIndex;
44  for (auto const &i : _myMap) {
45  auto mapping = i.second.get();
46  mapping->setIndex(ipar);
47  ipar += mapping->getNpar();
48  }
49  return ipar;
50 }
51 
52 void SimplePhotometryModel::offsetParams(Eigen::VectorXd const &delta) {
53  for (auto &i : _myMap) {
54  auto mapping = i.second.get();
55  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
56  }
57 }
58 
60  for (auto &i : _myMap) {
61  i.second->freezeErrorTransform();
62  }
63 }
64 
66  std::vector<unsigned> &indices) const {
67  auto mapping = findMapping(ccdImage);
68  if (indices.size() < mapping->getNpar()) indices.resize(mapping->getNpar());
69  indices[0] = mapping->getIndex();
70 }
71 
73  int total = 0;
74  for (auto &i : _myMap) {
75  total += i.second->getNpar();
76  }
77  return total;
78 }
79 
81  CcdImage const &ccdImage,
82  Eigen::VectorXd &derivatives) const {
83  auto mapping = findMapping(ccdImage);
84  mapping->computeParameterDerivatives(measuredStar, measuredStar.getInstFlux(), derivatives);
85 }
86 
88  for (auto &i : _myMap) {
89  stream << i.first << ": ";
90  i.second->dump(stream);
91  stream << ", ";
92  }
93 }
94 
96  auto i = _myMap.find(ccdImage.getHashKey());
97  if (i == _myMap.end())
99  "SimplePhotometryModel cannot find CcdImage " + ccdImage.getName());
100  return i->second.get();
101 }
102 
103 SimpleFluxModel::SimpleFluxModel(CcdImageList const &ccdImageList, double errorPedestal_)
104  : SimplePhotometryModel(ccdImageList, LOG_GET("jointcal.SimpleFluxModel"), errorPedestal_) {
105  for (auto const &ccdImage : ccdImageList) {
106  auto photoCalib = ccdImage->getPhotoCalib();
107  // Use the single-frame processing calibration from the PhotoCalib as the initial value.
108  auto transform = std::make_shared<FluxTransformSpatiallyInvariant>(photoCalib->getCalibrationMean());
109  _myMap.emplace(ccdImage->getHashKey(), std::make_unique<PhotometryMapping>(transform));
110  }
111  LOGLS_INFO(_log, "SimpleFluxModel got " << _myMap.size() << " ccdImage mappings.");
112 }
113 
114 double SimpleFluxModel::computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const {
115  return transform(ccdImage, measuredStar) - measuredStar.getFittedStar()->getFlux();
116 }
117 
118 double SimpleFluxModel::transform(CcdImage const &ccdImage, MeasuredStar const &star) const {
119  auto mapping = findMapping(ccdImage);
120  return mapping->transform(star, star.getInstFlux());
121 }
122 
124  auto mapping = findMapping(ccdImage);
125  double tempErr = tweakFluxError(star);
126  return mapping->transformError(star, star.getInstFlux(), tempErr);
127 }
128 
130  double calibration = (findMapping(ccdImage)->getParameters()[0]);
131  auto oldPhotoCalib = ccdImage.getPhotoCalib();
132  return std::make_unique<afw::image::PhotoCalib>(calibration, oldPhotoCalib->getCalibrationErr());
133 }
134 
135 SimpleMagnitudeModel::SimpleMagnitudeModel(CcdImageList const &ccdImageList, double errorPedestal_)
136  : SimplePhotometryModel(ccdImageList, LOG_GET("jointcal.SimpleMagnitudeModel"), errorPedestal_) {
137  for (auto const &ccdImage : ccdImageList) {
138  auto photoCalib = ccdImage->getPhotoCalib();
139  // Use the single-frame processing calibration from the PhotoCalib as the default.
140  double calib = magFromFlux(photoCalib->getCalibrationMean());
141  auto transform = std::make_shared<MagnitudeTransformSpatiallyInvariant>(calib);
142  _myMap.emplace(ccdImage->getHashKey(), std::make_unique<PhotometryMapping>(transform));
143  }
144  LOGLS_INFO(_log, "SimpleMagnitudeModel got " << _myMap.size() << " ccdImage mappings.");
145 }
146 
148  MeasuredStar const &measuredStar) const {
149  return transform(ccdImage, measuredStar) - measuredStar.getFittedStar()->getMag();
150 }
151 
153  auto mapping = findMapping(ccdImage);
154  return mapping->transform(star, star.getInstMag());
155 }
156 
158  auto mapping = findMapping(ccdImage);
159  double tempErr = tweakMagnitudeError(star);
160  return mapping->transformError(star, star.getInstMag(), tempErr);
161 }
162 
164  // NOTE: photocalib is defined as instFlux * calibration = flux
165  double calibration = instFluxFromMag(findMapping(ccdImage)->getParameters()[0]);
166  auto oldPhotoCalib = ccdImage.getPhotoCalib();
167  return std::make_unique<afw::image::PhotoCalib>(calibration, oldPhotoCalib->getCalibrationErr());
168 }
169 
170 } // namespace jointcal
171 } // namespace lsst
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
Key< int > calib
Relates transform(s) to their position in the fitting matrix and allows interaction with the transfor...
virtual Eigen::VectorXd getParameters()=0
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:79
int getTotalParameters() const override
Return the total number of parameters in this model.
SimpleMagnitudeModel(CcdImageList const &ccdImageList, double errorPedestal_=0)
std::shared_ptr< afw::image::PhotoCalib > toPhotoCalib(CcdImage const &ccdImage) const override
Return the mapping of ccdImage represented as a PhotoCalib.
T end(T... args)
void computeParameterDerivatives(MeasuredStar const &measuredStar, CcdImage const &ccdImage, Eigen::VectorXd &derivatives) const override
Compute the parametric derivatives of this model.
PhotometryMappingBase * findMapping(CcdImage const &ccdImage) const override
Return the mapping associated with this ccdImage.
T resize(T... args)
#define LOGLS_INFO(logger, message)
STL class.
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.
std::shared_ptr< afw::image::PhotoCalib > getPhotoCalib() const
Return the exposure&#39;s photometric calibration.
Definition: CcdImage.h:161
Class for a simple mapping implementing a generic AstrometryTransform.
SimpleFluxModel(CcdImageList const &ccdImageList, double errorPedestal_=0)
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
objects measured on actual images.
Definition: MeasuredStar.h:42
double transform(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux for measuredStar on ccdImage.
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 find(T... args)
T size(T... args)
#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).
double tweakMagnitudeError(jointcal::MeasuredStar const &measuredStar) const
Add a small magnitude offset to the "instrumental magnitude" error, in quadrature.
double computeResidual(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Compute the residual between the model applied to a star and its associated fittedStar.
T pow(T... args)
T emplace(T... args)
CcdImageKey getHashKey() const
Definition: CcdImage.h:152
void dump(std::ostream &stream=std::cout) const override
Dump the contents of the transforms, for debugging.
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
std::shared_ptr< FittedStar > getFittedStar() const
Definition: MeasuredStar.h:108
STL class.
#define LOG_GET(logger)
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
Photometric response model which has a single photometric factor per CcdImage.
double transformError(CcdImage const &ccdImage, MeasuredStar const &measuredStar) const override
Return the on-sky transformed flux uncertainty for measuredStar on ccdImage.
void freezeErrorTransform() override
Once this routine has been called, the error transform is not modified by offsetParams().