lsst.jointcal  16.0-17-g0bdc215
PhotometryFit.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <iomanip>
3 #include <algorithm>
4 #include <fstream>
5 
6 #include "Eigen/Sparse"
7 
8 #include "lsst/log/Log.h"
9 #include "lsst/pex/exceptions.h"
12 #include "lsst/jointcal/Chi2.h"
14 #include "lsst/jointcal/Gtransfo.h"
16 
17 namespace lsst {
18 namespace jointcal {
19 
20 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList,
21  Eigen::VectorXd &grad,
22  MeasuredStarList const *measuredStarList) const {
23  /**********************************************************************/
24  /* @note the math in this method and accumulateStatImageList() must be kept consistent,
25  * in terms of +/- convention, definition of model, etc. */
26  /**********************************************************************/
27 
28  /* this routine works in two different ways: either providing the
29  Ccd, of providing the MeasuredStarList. In the latter case, the
30  Ccd should match the one(s) in the list. */
31  if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
32 
33  unsigned nparModel = (_fittingModel) ? _photometryModel->getNpar(ccdImage) : 0;
34  unsigned nparFlux = (_fittingFluxes) ? 1 : 0;
35  unsigned nparTotal = nparModel + nparFlux;
36  std::vector<unsigned> indices(nparModel, -1);
37  if (_fittingModel) _photometryModel->getMappingIndices(ccdImage, indices);
38 
39  Eigen::VectorXd H(nparTotal); // derivative matrix
40  // current position in the Jacobian
41  unsigned kTriplets = tripletList.getNextFreeIndex();
42  const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
43 
44  for (auto const &measuredStar : catalog) {
45  if (!measuredStar->isValid()) continue;
46 // tweak the measurement errors
47 #ifdef FUTURE
48  TweakPhotomMeasurementErrors(inPos, *measuredStar, _fluxError);
49 #endif
50  H.setZero(); // we cannot be sure that all entries will be overwritten.
51 
52  double residual = _photometryModel->computeResidual(ccdImage, *measuredStar);
53  double inverseSigma = 1.0 / _photometryModel->transformError(ccdImage, *measuredStar);
54  double W = std::pow(inverseSigma, 2);
55 
56  if (_fittingModel) {
57  _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
58  for (unsigned k = 0; k < indices.size(); k++) {
59  unsigned l = indices[k];
60  tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
61  grad[l] += H[k] * W * residual;
62  }
63  }
64  if (_fittingFluxes) {
65  unsigned index = measuredStar->getFittedStar()->getIndexInMatrix();
66  // Note: H = dR/dFittedStarFlux == -1
67  tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
68  grad[index] += -1.0 * W * residual;
69  }
70  kTriplets += 1; // each measurement contributes 1 column in the Jacobian
71  }
72 
73  tripletList.setNextFreeIndex(kTriplets);
74 }
75 
76 void PhotometryFit::leastSquareDerivativesReference(FittedStarList const &fittedStarList,
77  TripletList &tripletList, Eigen::VectorXd &grad) const {
78  /**********************************************************************/
81  /**********************************************************************/
82 
83  // Derivatives of terms involving fitted and refstars only contribute if we are fitting fluxes.
84  if (!_fittingFluxes) return;
85  // Can't compute anything if there are no refStars.
86  if (_associations->refStarList.size() == 0) return;
87 
88  unsigned kTriplets = tripletList.getNextFreeIndex();
89 
90  for (auto const &fittedStar : fittedStarList) {
91  auto refStar = fittedStar->getRefStar();
92  if (refStar == nullptr) continue; // no contribution if no associated refstar
93  // TODO: Can we actually work with multiple filters at a time in this scheme?
94  // TODO: I feel we might only be able to fit one filter at a time, because these terms are
95  // independent of the ccdImage, so we don't have a specific filter.
96  // filter = ccdImage.getFilter();
97 
98  // W == inverseSigma^2
99 
100  double inverseSigma = 1.0 / _photometryModel->getRefError(*refStar);
101  // Residual is fittedStar - refStar for consistency with measurement terms.
102  double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
103 
104  unsigned index = fittedStar->getIndexInMatrix();
105  // Note: H = dR/dFittedStar == 1
106  tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
107  grad(index) += 1.0 * std::pow(inverseSigma, 2) * residual;
108  kTriplets += 1;
109  }
110  tripletList.setNextFreeIndex(kTriplets);
111 }
112 
113 void PhotometryFit::accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const {
114  /**********************************************************************/
117  /**********************************************************************/
118  for (auto const &ccdImage : ccdImageList) {
119  auto &catalog = ccdImage->getCatalogForFit();
120 
121  for (auto const &measuredStar : catalog) {
122  if (!measuredStar->isValid()) continue;
123  double sigma = _photometryModel->transformError(*ccdImage, *measuredStar);
124 #ifdef FUTURE
125  TweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
126 #endif
127  double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
128 
129  double chi2Val = std::pow(residual / sigma, 2);
130  accum.addEntry(chi2Val, 1, measuredStar);
131  } // end loop on measurements
132  }
133 }
134 
135 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum) const {
136  /**********************************************************************/
139  /**********************************************************************/
140 
141  FittedStarList &fittedStarList = _associations->fittedStarList;
142  for (auto const &fittedStar : fittedStarList) {
143  auto refStar = fittedStar->getRefStar();
144  if (refStar == nullptr) continue;
145  double sigma = _photometryModel->getRefError(*refStar);
146  double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
147  double chi2 = std::pow(residual / sigma, 2);
148  accum.addEntry(chi2, 1, fittedStar);
149  }
150 }
151 
153 
155 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar,
156  std::vector<unsigned> &indices) const {
157  indices.clear();
158  if (_fittingModel) {
159  _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
160  }
161  if (_fittingFluxes) {
162  std::shared_ptr<FittedStar const> const fs = measuredStar.getFittedStar();
163  unsigned fsIndex = fs->getIndexInMatrix();
164  indices.push_back(fsIndex);
165  }
166 }
167 
169  _whatToFit = whatToFit;
170  LOGLS_INFO(_log, "assignIndices: now fitting: " << whatToFit);
171  _fittingModel = (_whatToFit.find("Model") != std::string::npos);
172  _fittingFluxes = (_whatToFit.find("Fluxes") != std::string::npos);
173  // When entering here, we assume that whatToFit has already been interpreted.
174 
175  _nParModel = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
176  unsigned ipar = _nParModel;
177 
178  if (_fittingFluxes) {
179  for (auto &fittedStar : _associations->fittedStarList) {
180  // the parameter layout here is used also
181  // - when filling the derivatives
182  // - when updating (offsetParams())
183  // - in getIndicesOfMeasuredStar
184  fittedStar->setIndexInMatrix(ipar);
185  ipar += 1;
186  }
187  }
188  _nParFluxes = ipar - _nParModel;
189  _nParTot = ipar;
191  "nParameters total: " << _nParTot << " model: " << _nParModel << " fluxes: " << _nParFluxes);
192 }
193 
194 void PhotometryFit::offsetParams(Eigen::VectorXd const &delta) {
195  if (delta.size() != _nParTot)
197  "PhotometryFit::offsetParams : the provided vector length is not compatible with "
198  "the current whatToFit setting");
199  if (_fittingModel) _photometryModel->offsetParams(delta);
200 
201  if (_fittingFluxes) {
202  for (auto &fittedStar : _associations->fittedStarList) {
203  // the parameter layout here is used also
204  // - when filling the derivatives
205  // - when assigning indices (assignIndices())
206  unsigned index = fittedStar->getIndexInMatrix();
207  _photometryModel->offsetFittedStar(*fittedStar, delta(index));
208  }
209  }
210 }
211 
213  std::ofstream ofile(baseName.c_str());
214  std::string separator = "\t";
215 
216  ofile << "#id" << separator << "xccd" << separator << "yccd" << separator;
217  ofile << "mag" << separator << "instMag" << separator << "instMagErr" << separator;
218  ofile << "instFlux" << separator << "instFluxErr" << separator;
219  ofile << "inputFlux" << separator << "inputFluxErr" << separator;
220  ofile << "transformedFlux" << separator << "transformedFluxErr" << separator;
221  ofile << "fittedFlux" << separator;
222  ofile << "mjd" << separator << "color" << separator;
223  ofile << "fsindex" << separator;
224  ofile << "ra" << separator << "dec" << separator;
225  ofile << "chi2" << separator << "nm" << separator;
226  ofile << "chip" << separator << "visit" << separator << std::endl;
227 
228  ofile << "#id in source catalog" << separator << "coordinates in CCD" << separator << separator;
229  ofile << "fitted magnitude" << separator << "measured magnitude" << separator
230  << "measured magnitude error" << separator;
231  ofile << "measured instrument flux (ADU)" << separator << "measured instrument flux error" << separator;
232  ofile << "measured flux (maggies)" << separator << "measured flux error" << separator;
233  ofile << separator << separator;
234  ofile << "fitted flux (maggies)" << separator;
235  ofile << "modified Julian date of the measurement" << separator << "currently unused" << separator;
236  ofile << "unique index of the fittedStar" << separator;
237  ofile << "on-sky position of fitted star" << separator << separator;
238  ofile << "contribution to Chi2 (1 dof)" << separator << "number of measurements of this FittedStar"
239  << separator;
240  ofile << "chip id" << separator << "visit id" << std::endl;
241 
242  const CcdImageList &ccdImageList = _associations->getCcdImageList();
243  for (auto const &ccdImage : ccdImageList) {
244  const MeasuredStarList &cat = ccdImage->getCatalogForFit();
245  for (auto const &measuredStar : cat) {
246  if (!measuredStar->isValid()) continue;
247  double sigma = _photometryModel->transformError(*ccdImage, *measuredStar);
248 #ifdef FUTURE
249  tweakPhotomMeasurementErrors(inPos, measuredStar, _fluxError);
250 #endif
251  double flux = _photometryModel->transform(*ccdImage, *measuredStar);
252  double fluxErr = _photometryModel->transformError(*ccdImage, *measuredStar);
253  double jd = ccdImage->getMjd();
254  std::shared_ptr<FittedStar const> const fittedStar = measuredStar->getFittedStar();
255  double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
256  double chi2Val = std::pow(residual / sigma, 2);
257 
258  ofile << std::setprecision(9);
259  ofile << measuredStar->getId() << separator << measuredStar->x << separator << measuredStar->y
260  << separator;
261  ofile << fittedStar->getMag() << separator << measuredStar->getInstMag() << separator
262  << measuredStar->getInstMagErr() << separator;
263  ofile << measuredStar->getInstFlux() << separator << measuredStar->getInstFluxErr() << separator;
264  ofile << measuredStar->getFlux() << separator << measuredStar->getFluxErr() << separator;
265  ofile << flux << separator << fluxErr << separator << fittedStar->getFlux() << separator;
266  ofile << jd << separator << fittedStar->color << separator;
267  ofile << fittedStar->getIndexInMatrix() << separator;
268  ofile << fittedStar->x << separator << fittedStar->y << separator;
269  ofile << chi2Val << separator << fittedStar->getMeasurementCount() << separator;
270  ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
271  } // loop on measurements in image
272  } // loop on images
273 }
274 
276  std::ofstream ofile(baseName.c_str());
277  std::string separator = "\t";
278 
279  ofile << "#ra" << separator << "dec " << separator;
280  ofile << "mag" << separator << "color" << separator;
281  ofile << "refFlux" << separator << "refFluxErr" << separator;
282  ofile << "fittedFlux" << separator << "fittedFluxErr" << separator;
283  ofile << "fsindex" << separator << "chi2" << separator << "nm" << std::endl;
284 
285  ofile << "#coordinates of fittedStar" << separator << separator;
286  ofile << "magnitude" << separator << "currently unused" << separator;
287  ofile << "refStar flux (maggies)" << separator << "refStar fluxErr" << separator;
288  ofile << "fittedStar flux (maggies)" << separator << "fittedStar fluxErr" << separator;
289  ofile << "unique index of the fittedStar" << separator << "refStar contribution to Chi2 (1 dof)"
290  << separator << "number of measurements of this FittedStar" << std::endl;
291 
292  // The following loop is heavily inspired from PhotometryFit::computeChi2()
293  const FittedStarList &fittedStarList = _associations->fittedStarList;
294  for (auto const &fittedStar : fittedStarList) {
295  const RefStar *refStar = fittedStar->getRefStar();
296  if (refStar == nullptr) continue;
297 
298  double chi2 = std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
299 
300  ofile << std::setprecision(9);
301  ofile << fittedStar->x << separator << fittedStar->y << separator;
302  ofile << fittedStar->getMag() << separator << fittedStar->color << separator;
303  ofile << refStar->getFlux() << separator << refStar->getFluxErr() << separator;
304  ofile << fittedStar->getFlux() << separator << fittedStar->getFluxErr() << separator;
305  ofile << fittedStar->getIndexInMatrix() << separator << chi2 << separator
306  << fittedStar->getMeasurementCount() << std::endl;
307  } // loop on FittedStars
308 }
309 
310 } // namespace jointcal
311 } // namespace lsst
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:16
double getFluxErr(size_t filter) const
reference fluxErr in a given filter
Definition: RefStar.h:44
void saveChi2RefContributions(std::string const &baseName) const override
Save a CSV file containing residuals of reference terms.
T endl(T... args)
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:118
afw::table::Key< double > sigma
#define LOGLS_INFO(logger, message)
STL class.
T push_back(T... args)
STL class.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
Class for a simple mapping implementing a generic Gtransfo.
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:99
T clear(T... args)
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:136
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
T find(T... args)
#define LSST_EXCEPT(type,...)
T pow(T... args)
T c_str(T... args)
#define LOGLS_DEBUG(logger, message)
double getFlux(size_t filter) const
reference flux in a given filter
Definition: RefStar.h:42
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:23
T setprecision(T... args)
void saveChi2MeasContributions(std::string const &baseName) const override
Save a CSV file containing residuals of measurement terms.