LSST Applications  21.0.0-147-g0e635eb1+1acddb5be5,22.0.0+052faf71bd,22.0.0+1ea9a8b2b2,22.0.0+6312710a6c,22.0.0+729191ecac,22.0.0+7589c3a021,22.0.0+9f079a9461,22.0.1-1-g7d6de66+b8044ec9de,22.0.1-1-g87000a6+536b1ee016,22.0.1-1-g8e32f31+6312710a6c,22.0.1-10-gd060f87+016f7cdc03,22.0.1-12-g9c3108e+df145f6f68,22.0.1-16-g314fa6d+c825727ab8,22.0.1-19-g93a5c75+d23f2fb6d8,22.0.1-19-gb93eaa13+aab3ef7709,22.0.1-2-g8ef0a89+b8044ec9de,22.0.1-2-g92698f7+9f079a9461,22.0.1-2-ga9b0f51+052faf71bd,22.0.1-2-gac51dbf+052faf71bd,22.0.1-2-gb66926d+6312710a6c,22.0.1-2-gcb770ba+09e3807989,22.0.1-20-g32debb5+b8044ec9de,22.0.1-23-gc2439a9a+fb0756638e,22.0.1-3-g496fd5d+09117f784f,22.0.1-3-g59f966b+1e6ba2c031,22.0.1-3-g849a1b8+f8b568069f,22.0.1-3-gaaec9c0+c5c846a8b1,22.0.1-32-g5ddfab5d3+60ce4897b0,22.0.1-4-g037fbe1+64e601228d,22.0.1-4-g8623105+b8044ec9de,22.0.1-5-g096abc9+d18c45d440,22.0.1-5-g15c806e+57f5c03693,22.0.1-7-gba73697+57f5c03693,master-g6e05de7fdc+c1283a92b8,master-g72cdda8301+729191ecac,w.2021.39
LSST Data Management Base Package
PhotometryFit.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 <iomanip>
27 #include <algorithm>
28 #include <fstream>
29 #include <cmath>
30 
31 #include "Eigen/Sparse"
32 
33 #include "lsst/log/Log.h"
34 #include "lsst/pex/exceptions.h"
37 #include "lsst/jointcal/Chi2.h"
41 
42 namespace lsst {
43 namespace jointcal {
44 
45 void PhotometryFit::leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList,
46  Eigen::VectorXd &grad,
47  MeasuredStarList const *measuredStarList) const {
48  /**********************************************************************/
49  /* @note the math in this method and accumulateStatImageList() must be kept consistent,
50  * in terms of +/- convention, definition of model, etc. */
51  /**********************************************************************/
52 
53  /* this routine works in two different ways: either providing the
54  Ccd, of providing the MeasuredStarList. In the latter case, the
55  Ccd should match the one(s) in the list. */
56  if (measuredStarList) assert(&(measuredStarList->front()->getCcdImage()) == &ccdImage);
57 
58  std::size_t nparModel = (_fittingModel) ? _photometryModel->getNpar(ccdImage) : 0;
59  std::size_t nparFlux = (_fittingFluxes) ? 1 : 0;
60  std::size_t nparTotal = nparModel + nparFlux;
61  IndexVector indices(nparModel, -1);
62  if (_fittingModel) _photometryModel->getMappingIndices(ccdImage, indices);
63 
64  Eigen::VectorXd H(nparTotal); // derivative matrix
65  // current position in the Jacobian
66  Eigen::Index kTriplets = tripletList.getNextFreeIndex();
67  const MeasuredStarList &catalog = (measuredStarList) ? *measuredStarList : ccdImage.getCatalogForFit();
68 
69  for (auto const &measuredStar : catalog) {
70  if (!measuredStar->isValid()) continue;
71  H.setZero(); // we cannot be sure that all entries will be overwritten.
72 
73  double residual = _photometryModel->computeResidual(ccdImage, *measuredStar);
74  double inverseSigma = 1.0 / _photometryModel->transformError(ccdImage, *measuredStar);
75  double W = std::pow(inverseSigma, 2);
76 
77  if (_fittingModel) {
78  _photometryModel->computeParameterDerivatives(*measuredStar, ccdImage, H);
79  for (std::size_t k = 0; k < indices.size(); k++) {
80  Eigen::Index l = indices[k];
81  tripletList.addTriplet(l, kTriplets, H[k] * inverseSigma);
82  grad[l] += H[k] * W * residual;
83  }
84  }
85  if (_fittingFluxes) {
86  Eigen::Index index = measuredStar->getFittedStar()->getIndexInMatrix();
87  // Note: H = dR/dFittedStarFlux == -1
88  tripletList.addTriplet(index, kTriplets, -1.0 * inverseSigma);
89  grad[index] += -1.0 * W * residual;
90  }
91  kTriplets += 1; // each measurement contributes 1 column in the Jacobian
92  }
93 
94  tripletList.setNextFreeIndex(kTriplets);
95 }
96 
97 void PhotometryFit::leastSquareDerivativesReference(FittedStarList const &fittedStarList,
98  TripletList &tripletList, Eigen::VectorXd &grad) const {
99  /**********************************************************************/
102  /**********************************************************************/
103 
104  // Derivatives of terms involving fitted and refstars only contribute if we are fitting fluxes.
105  if (!_fittingFluxes) return;
106  // Can't compute anything if there are no refStars.
107  if (_associations->refStarList.size() == 0) return;
108 
109  Eigen::Index kTriplets = tripletList.getNextFreeIndex();
110 
111  for (auto const &fittedStar : fittedStarList) {
112  auto refStar = fittedStar->getRefStar();
113  if (refStar == nullptr) continue; // no contribution if no associated refstar
114  // TODO: Can we actually work with multiple filters at a time in this scheme?
115  // TODO: I feel we might only be able to fit one filter at a time, because these terms are
116  // independent of the ccdImage, so we don't have a specific filter.
117  // filter = ccdImage.getFilter();
118 
119  // W == inverseSigma^2
120 
121  double inverseSigma = 1.0 / _photometryModel->getRefError(*refStar);
122  // Residual is fittedStar - refStar for consistency with measurement terms.
123  double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
124 
125  Eigen::Index index = fittedStar->getIndexInMatrix();
126  // Note: H = dR/dFittedStar == 1
127  tripletList.addTriplet(index, kTriplets, 1.0 * inverseSigma);
128  grad(index) += 1.0 * std::pow(inverseSigma, 2) * residual;
129  kTriplets += 1;
130  }
131  tripletList.setNextFreeIndex(kTriplets);
132 }
133 
134 void PhotometryFit::accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const {
135  /**********************************************************************/
138  /**********************************************************************/
139  for (auto const &ccdImage : ccdImageList) {
140  auto &catalog = ccdImage->getCatalogForFit();
141 
142  for (auto const &measuredStar : catalog) {
143  if (!measuredStar->isValid()) continue;
144  double sigma = _photometryModel->transformError(*ccdImage, *measuredStar);
145  double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
146 
147  double chi2Val = std::pow(residual / sigma, 2);
148  accum.addEntry(chi2Val, 1, measuredStar);
149  } // end loop on measurements
150  }
151 }
152 
153 void PhotometryFit::accumulateStatRefStars(Chi2Accumulator &accum) const {
154  /**********************************************************************/
157  /**********************************************************************/
158 
159  FittedStarList &fittedStarList = _associations->fittedStarList;
160  for (auto const &fittedStar : fittedStarList) {
161  auto refStar = fittedStar->getRefStar();
162  if (refStar == nullptr) continue;
163  double sigma = _photometryModel->getRefError(*refStar);
164  double residual = _photometryModel->computeRefResidual(*fittedStar, *refStar);
165  double chi2 = std::pow(residual / sigma, 2);
166  accum.addEntry(chi2, 1, fittedStar);
167  }
168 }
169 
171 
173 void PhotometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const {
174  indices.clear();
175  if (_fittingModel) {
176  _photometryModel->getMappingIndices(measuredStar.getCcdImage(), indices);
177  }
178  if (_fittingFluxes) {
179  std::shared_ptr<FittedStar const> const fs = measuredStar.getFittedStar();
180  Eigen::Index fsIndex = fs->getIndexInMatrix();
181  indices.push_back(fsIndex);
182  }
183 }
184 
186  _whatToFit = whatToFit;
187  LOGLS_INFO(_log, "assignIndices: now fitting: " << whatToFit);
188  _fittingModel = (_whatToFit.find("Model") != std::string::npos);
189  _fittingFluxes = (_whatToFit.find("Fluxes") != std::string::npos);
190  // When entering here, we assume that whatToFit has already been interpreted.
191 
192  _nModelParams = (_fittingModel) ? _photometryModel->assignIndices(whatToFit, 0) : 0;
193  std::size_t ipar = _nModelParams;
194 
195  if (_fittingFluxes) {
196  for (auto &fittedStar : _associations->fittedStarList) {
197  // the parameter layout here is used also
198  // - when filling the derivatives
199  // - when updating (offsetParams())
200  // - in getIndicesOfMeasuredStar
201  fittedStar->setIndexInMatrix(ipar);
202  ipar += 1;
203  }
204  }
205  _nStarParams = ipar - _nModelParams;
206  _nTotal = ipar;
207  LOGLS_DEBUG(_log, "nParameters total: " << _nTotal << " model: " << _nModelParams
208  << " fluxes: " << _nStarParams);
209 }
210 
211 void PhotometryFit::offsetParams(Eigen::VectorXd const &delta) {
212  if (delta.size() != _nTotal)
214  "PhotometryFit::offsetParams : the provided vector length is not compatible with "
215  "the current whatToFit setting");
216  if (_fittingModel) _photometryModel->offsetParams(delta);
217 
218  if (_fittingFluxes) {
219  for (auto &fittedStar : _associations->fittedStarList) {
220  // the parameter layout here is used also
221  // - when filling the derivatives
222  // - when assigning indices (assignIndices())
223  Eigen::Index index = fittedStar->getIndexInMatrix();
224  _photometryModel->offsetFittedStar(*fittedStar, delta(index));
225  }
226  }
227 }
228 
230  std::ofstream ofile(filename.c_str());
231  std::string separator = "\t";
232 
233  ofile << "id" << separator << "xccd" << separator << "yccd" << separator;
234  ofile << "mag" << separator << "instMag" << separator << "instMagErr" << separator;
235  ofile << "instFlux" << separator << "instFluxErr" << separator;
236  ofile << "inputFlux" << separator << "inputFluxErr" << separator;
237  ofile << "transformedFlux" << separator << "transformedFluxErr" << separator;
238  ofile << "fittedFlux" << separator;
239  ofile << "epoch" << separator;
240  ofile << "fsindex" << separator;
241  ofile << "ra" << separator << "dec" << separator;
242  ofile << "chi2" << separator << "nm" << separator;
243  ofile << "detector" << separator << "visit" << separator << std::endl;
244 
245  ofile << "id in source catalog" << separator << "coordinates in CCD" << separator << separator;
246  ofile << "fitted magnitude" << separator << "measured magnitude" << separator
247  << "measured magnitude error" << separator;
248  ofile << "measured instrumental flux (ADU)" << separator << "measured instrument flux error" << separator;
249  ofile << "measured flux (nJy)" << separator << "measured flux error" << separator;
250  ofile << separator << separator;
251  ofile << "fitted flux (nJy)" << separator;
252  ofile << "Julian Epoch Year of the measurement" << separator;
253  ofile << "unique index of the fittedStar" << separator;
254  ofile << "on-sky position of fitted star" << separator << separator;
255  ofile << "contribution to chi2 (1 dof)" << separator << "number of measurements of this FittedStar"
256  << separator;
257  ofile << "detector id" << separator << "visit id" << std::endl;
258 
259  const CcdImageList &ccdImageList = _associations->getCcdImageList();
260  for (auto const &ccdImage : ccdImageList) {
261  const MeasuredStarList &cat = ccdImage->getCatalogForFit();
262  for (auto const &measuredStar : cat) {
263  if (!measuredStar->isValid()) continue;
264 
265  double instFluxErr = _photometryModel->tweakFluxError(*measuredStar);
266  double flux = _photometryModel->transform(*ccdImage, *measuredStar);
267  double fluxErr = _photometryModel->transformError(*ccdImage, *measuredStar);
268  std::shared_ptr<FittedStar const> const fittedStar = measuredStar->getFittedStar();
269  double residual = _photometryModel->computeResidual(*ccdImage, *measuredStar);
270  double chi2Val = std::pow(residual / fluxErr, 2);
271 
272  ofile << std::setprecision(9);
273  ofile << measuredStar->getId() << separator << measuredStar->x << separator << measuredStar->y
274  << separator;
275  ofile << fittedStar->getMag() << separator << measuredStar->getInstMag() << separator
276  << measuredStar->getInstMagErr() << separator;
277  ofile << measuredStar->getInstFlux() << separator << instFluxErr << separator;
278  ofile << measuredStar->getFlux() << separator << measuredStar->getFluxErr() << separator;
279  ofile << flux << separator << fluxErr << separator << fittedStar->getFlux() << separator;
280  ofile << ccdImage->getEpoch() << separator;
281  ofile << fittedStar->getIndexInMatrix() << separator;
282  ofile << fittedStar->x << separator << fittedStar->y << separator;
283  ofile << chi2Val << separator << fittedStar->getMeasurementCount() << separator;
284  ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
285  } // loop on measurements in image
286  } // loop on images
287 }
288 
290  std::ofstream ofile(filename.c_str());
291  std::string separator = "\t";
292 
293  ofile << "ra" << separator << "dec " << separator;
294  ofile << "mag" << separator;
295  ofile << "refFlux" << separator << "refFluxErr" << separator;
296  ofile << "fittedFlux" << separator << "fittedFluxErr" << separator;
297  ofile << "fsindex" << separator << "chi2" << separator << "nm" << std::endl;
298 
299  ofile << "coordinates of fittedStar" << separator << separator;
300  ofile << "magnitude" << separator;
301  ofile << "refStar flux (nJy)" << separator << "refStar fluxErr" << separator;
302  ofile << "fittedStar flux (nJy)" << separator << "fittedStar fluxErr" << separator;
303  ofile << "unique index of the fittedStar" << separator << "refStar contribution to chi2 (1 dof)"
304  << separator << "number of measurements of this FittedStar" << std::endl;
305 
306  // The following loop is heavily inspired from PhotometryFit::computeChi2()
307  const FittedStarList &fittedStarList = _associations->fittedStarList;
308  for (auto const &fittedStar : fittedStarList) {
309  const RefStar *refStar = fittedStar->getRefStar();
310  if (refStar == nullptr) continue;
311 
312  double chi2 = std::pow(((fittedStar->getFlux() - refStar->getFlux()) / refStar->getFluxErr()), 2);
313 
314  ofile << std::setprecision(9);
315  ofile << fittedStar->x << separator << fittedStar->y << separator;
316  ofile << fittedStar->getMag() << separator;
317  ofile << refStar->getFlux() << separator << refStar->getFluxErr() << separator;
318  ofile << fittedStar->getFlux() << separator << fittedStar->getFluxErr() << separator;
319  ofile << fittedStar->getIndexInMatrix() << separator << chi2 << separator
320  << fittedStar->getMeasurementCount() << std::endl;
321  } // loop on FittedStars
322 }
323 
324 } // namespace jointcal
325 } // namespace lsst
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
afw::table::Key< double > sigma
Definition: GaussianPsf.cc:49
LSST DM logging module built on log4cxx.
#define LOGLS_INFO(logger, message)
Log a info-level message using an iostream-based interface.
Definition: Log.h:639
#define LOGLS_DEBUG(logger, message)
Log a debug-level message using an iostream-based interface.
Definition: Log.h:619
T c_str(T... args)
double getFlux() const
Definition: BaseStar.h:98
double getFluxErr() const
Definition: BaseStar.h:102
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:116
Eigen::Index _nModelParams
Definition: FitterBase.h:168
std::shared_ptr< Associations > _associations
Definition: FitterBase.h:163
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:151
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
void saveChi2RefContributions(std::string const &filename) const override
Save a CSV file containing residuals of reference terms.
void saveChi2MeasContributions(std::string const &filename) const override
Save a CSV file containing residuals of measurement terms.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
Objects used as position/flux anchors (e.g.
Definition: RefStar.h:45
Reports invalid arguments.
Definition: Runtime.h:66
T clear(T... args)
T endl(T... args)
T find(T... args)
std::list< std::shared_ptr< CcdImage > > CcdImageList
Definition: CcdImage.h:46
A base class for image defects.
T pow(T... args)
T push_back(T... args)
T setprecision(T... args)