LSST Applications g180d380827+770a9040cc,g2079a07aa2+86d27d4dc4,g2305ad1205+09cfdadad9,g2bbee38e9b+c6a8a0fb72,g337abbeb29+c6a8a0fb72,g33d1c0ed96+c6a8a0fb72,g3a166c0a6a+c6a8a0fb72,g3ddfee87b4+1ea5e09c42,g48712c4677+7e2ea9cd42,g487adcacf7+301d09421d,g50ff169b8f+96c6868917,g52b1c1532d+585e252eca,g591dd9f2cf+96fcb956a6,g64a986408d+23540ee355,g858d7b2824+23540ee355,g864b0138d7+aa38e45daa,g95921f966b+d83dc58ecd,g991b906543+23540ee355,g99cad8db69+7f13b58a93,g9c22b2923f+e2510deafe,g9ddcbc5298+9a081db1e4,ga1e77700b3+03d07e1c1f,gb0e22166c9+60f28cb32d,gb23b769143+23540ee355,gba4ed39666+c2a2e4ac27,gbb8dafda3b+49e7449578,gbd998247f1+585e252eca,gc120e1dc64+1bbfa184e1,gc28159a63d+c6a8a0fb72,gc3e9b769f7+385ea95214,gcf0d15dbbd+1ea5e09c42,gdaeeff99f8+f9a426f77a,ge6526c86ff+1bccc98490,ge79ae78c31+c6a8a0fb72,gee10cc3b42+585e252eca,w.2024.18
LSST Data Management Base Package
Loading...
Searching...
No Matches
AstrometryFit.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 <algorithm>
26#include <fstream>
27#include <iomanip>
28#include <iostream>
29#include <utility>
30
31#include "Eigen/Sparse"
32
33#include "lsst/log/Log.h"
34#include "lsst/pex/exceptions.h"
37#include "lsst/jointcal/Chi2.h"
43
44namespace {
54double computeProjectedRefStarChi2(const lsst::jointcal::FatPoint& refStar) {
55 double det = refStar.vx * refStar.vy - std::pow(refStar.vxy, 2);
56 double wxx = refStar.vy / det;
57 double wyy = refStar.vx / det;
58 double wxy = -refStar.vxy / det;
59 return wxx * std::pow(refStar.x, 2) + 2 * wxy * refStar.x * refStar.y + wyy * std::pow(refStar.y, 2);
60}
61} // namespace
62
63namespace lsst {
64namespace jointcal {
65
67 std::shared_ptr<AstrometryModel> astrometryModel, double posError)
68 : FitterBase(std::move(associations)),
69 _astrometryModel(std::move(astrometryModel)),
70 _epoch(_associations->getEpoch()),
71 _posError(posError) {
72 _log = LOG_GET("lsst.jointcal.AstrometryFit");
73}
74
75/* ! this routine is used in 3 instances: when computing
76the derivatives, when computing the Chi2, when filling a tuple.
77*/
78Point AstrometryFit::transformFittedStar(FittedStar const &fittedStar, AstrometryTransform const &sky2TP,
79 double deltaYears) const {
80 Point fittedStarInTP;
81 if (fittedStar.getRefStar()) {
82 // PM data is on-sky, not in the tangent plane
83 Point temp = fittedStar.getRefStar()->applyProperMotion(fittedStar, deltaYears);
84 fittedStarInTP = sky2TP.apply(temp);
85 } else {
86 fittedStarInTP = sky2TP.apply(fittedStar);
87 }
88 return fittedStarInTP;
89}
90
94static void tweakAstromMeasurementErrors(FatPoint &P, MeasuredStar const &Ms, double error) {
95 static bool called = false;
96 static double increment = 0;
97 if (!called) {
98 increment = std::pow(error, 2); // was in Preferences
99 called = true;
100 }
101 P.vx += increment;
102 P.vy += increment;
103}
104
105// we could consider computing the chi2 here.
106// (although it is not extremely useful)
108 Eigen::VectorXd &fullGrad,
109 MeasuredStarList const *msList) const {
110 /**********************************************************************/
111 /* @note the math in this method and accumulateStatImage() must be kept consistent,
112 * in terms of +/- convention, definition of model, etc. */
113 /**********************************************************************/
114
115 /* Setup */
116 /* this routine works in two different ways: either providing the
117 ccdImage, of providing the MeasuredStarList. In the latter case, the
118 ccdImage should match the one(s) in the list. */
119 if (msList) assert(&(msList->front()->getCcdImage()) == &ccdImage);
120
121 // get the Mapping
122 const AstrometryMapping *mapping = _astrometryModel->getMapping(ccdImage);
123 // count parameters
124 std::size_t npar_mapping = (_fittingDistortions) ? mapping->getNpar() : 0;
125 std::size_t npar_pos = (_fittingPos) ? 2 : 0;
126 std::size_t npar_pm = (_fittingPM) ? 2 : 0;
127 std::size_t npar_tot = npar_mapping + npar_pos + npar_pm;
128 // if (npar_tot == 0) this CcdImage does not contribute
129 // any constraint to the fit, so :
130 if (npar_tot == 0) return;
131 IndexVector indices(npar_tot, -1);
132 if (_fittingDistortions) mapping->getMappingIndices(indices);
133
134 // FittedStar is "observed" epoch, MeasuredStar is "baseline"
135 double deltaYears = _epoch - ccdImage.getEpoch();
136 // transformation from sky to TP
137 auto sky2TP = _astrometryModel->getSkyToTangentPlane(ccdImage);
138 // reserve matrices once for all measurements
140 // the shape of H (et al) is required this way in order to be able to
141 // separate derivatives along x and y as vectors.
142 Eigen::MatrixX2d H(npar_tot, 2), halpha(npar_tot, 2), HW(npar_tot, 2);
143 Eigen::Matrix2d transW(2, 2);
144 Eigen::Matrix2d alpha(2, 2);
145 Eigen::VectorXd grad(npar_tot);
146 // current position in the Jacobian
147 Eigen::Index kTriplets = tripletList.getNextFreeIndex();
148 const MeasuredStarList &catalog = (msList) ? *msList : ccdImage.getCatalogForFit();
149
150 for (auto &i : catalog) {
151 const MeasuredStar &ms = *i;
152 if (!ms.isValid()) continue;
153 // tweak the measurement errors
154 FatPoint inPos = ms;
155 tweakAstromMeasurementErrors(inPos, ms, _posError);
156 H.setZero(); // we cannot be sure that all entries will be overwritten.
157 FatPoint outPos;
158 // should *not* fill H if whatToFit excludes mapping parameters.
159 if (_fittingDistortions)
160 mapping->computeTransformAndDerivatives(inPos, outPos, H);
161 else
162 mapping->transformPosAndErrors(inPos, outPos);
163
164 std::size_t ipar = npar_mapping;
165 double det = outPos.vx * outPos.vy - std::pow(outPos.vxy, 2);
166 if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
167 LOGLS_WARN(_log, "Inconsistent measurement errors: drop measurement at "
168 << Point(ms) << " in image " << ccdImage.getName());
169 continue;
170 }
171 transW(0, 0) = outPos.vy / det;
172 transW(1, 1) = outPos.vx / det;
173 transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
174 // compute alpha, a triangular square root
175 // of transW (i.e. a Cholesky factor)
176 alpha(0, 0) = sqrt(transW(0, 0));
177 // checked that alpha*alphaT = transW
178 alpha(1, 0) = transW(0, 1) / alpha(0, 0);
179 // DB - I think that the next line is equivalent to : alpha(1,1) = 1./sqrt(outPos.vy)
180 // PA - seems correct !
181 alpha(1, 1) = 1. / sqrt(det * transW(0, 0));
182 alpha(0, 1) = 0;
183
185
186 Point fittedStarInTP = transformFittedStar(*fs, *sky2TP, deltaYears);
187
188 // compute derivative of TP position w.r.t sky position ....
189 if (npar_pos > 0) // ... if actually fitting FittedStar position
190 {
191 sky2TP->computeDerivative(*fs, dypdy, 1e-3);
192 // sign checked
193 // TODO Still have to check with non trivial non-diagonal terms
194 H(npar_mapping, 0) = -dypdy.A11();
195 H(npar_mapping + 1, 0) = -dypdy.A12();
196 H(npar_mapping, 1) = -dypdy.A21();
197 H(npar_mapping + 1, 1) = -dypdy.A22();
198 indices[npar_mapping] = fs->getIndexInMatrix();
199 indices.at(npar_mapping + 1) = fs->getIndexInMatrix() + 1;
200 ipar += npar_pos;
201 }
202 /* only consider proper motions of objects allowed to move,
203 unless the fit is going to be degenerate */
204 // TODO: left as reference for when we implement PM fitting
205 // TODO: mjd here would become either deltaYears or maybe associations.epoch? Check the math first!
206 // if (_fittingPM && fs->mightMove) {
207 // H(ipar, 0) = -mjd; // Sign unchecked but consistent with above
208 // H(ipar + 1, 1) = -mjd;
209 // indices[ipar] = fs->getIndexInMatrix() + 2;
210 // indices[ipar + 1] = fs->getIndexInMatrix() + 3;
211 // ipar += npar_pm;
212 // }
213
214 // We can now compute the residual
215 Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
216
217 // do not write grad = H*transW*res to avoid
218 // dynamic allocation of a temporary
219 halpha = H * alpha;
220 HW = H * transW;
221 grad = HW * res;
222 // now feed in triplets and fullGrad
223 for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
224 for (std::size_t ic = 0; ic < 2; ++ic) {
225 double val = halpha(ipar, ic);
226 if (val == 0) continue;
227 tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
228 }
229 fullGrad(indices[ipar]) += grad(ipar);
230 }
231 kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
232 } // end loop on measurements
233 tripletList.setNextFreeIndex(kTriplets);
234}
235
237 TripletList &tripletList,
238 Eigen::VectorXd &fullGrad) const {
239 /**********************************************************************/
240 /* @note the math in this method and accumulateStatRefStars() must be kept consistent,
241 * in terms of +/- convention, definition of model, etc. */
242 /**********************************************************************/
243
244 /* We compute here the derivatives of the terms involving fitted
245 stars and reference stars. They only provide contributions if we
246 are fitting positions: */
247 if (!_fittingPos) return;
248 /* the other case where the accumulation of derivatives stops
249 here is when there are no RefStars */
250 if (_associations->refStarList.size() == 0) return;
251 Eigen::Matrix2d W(2, 2);
252 Eigen::Matrix2d alpha(2, 2);
253 Eigen::Matrix2d H(2, 2), halpha(2, 2), HW(2, 2);
255 Eigen::Vector2d res, grad;
256 Eigen::Index indices[2];
257 Eigen::Index kTriplets = tripletList.getNextFreeIndex();
258 /* We cannot use the spherical coordinates directly to evaluate
259 Euclidean distances, we have to use a projector on some plane in
260 order to express least squares. Not projecting could lead to a
261 disaster around the poles or across alpha=0. So we need a
262 projector. We construct a projector and will change its
263 projection point at every object */
265 for (auto const &i : fittedStarList) {
266 const FittedStar &fs = *i;
267 auto rs = fs.getRefStar();
268 if (rs == nullptr) continue;
269 proj.setTangentPoint(fs);
270 // fs projects to (0,0), no need to compute its transform.
271 FatPoint rsProj;
272 proj.transformPosAndErrors(*rs, rsProj);
273 // Compute the derivative of the projector to incorporate its effects on the errors.
274 proj.computeDerivative(fs, der, 1e-4);
275 // sign checked. TODO check that the off-diagonal terms are OK.
276 H(0, 0) = -der.A11();
277 H(1, 0) = -der.A12();
278 H(0, 1) = -der.A21();
279 H(1, 1) = -der.A22();
280 // TO DO : account for proper motions.
281 double det = rsProj.vx * rsProj.vy - std::pow(rsProj.vxy, 2);
282 if (rsProj.vx <= 0 || rsProj.vy <= 0 || det <= 0) {
283 LOGLS_WARN(_log, "RefStar error matrix not positive definite for: " << *rs);
284 continue;
285 }
286 W(0, 0) = rsProj.vy / det;
287 W(0, 1) = W(1, 0) = -rsProj.vxy / det;
288 W(1, 1) = rsProj.vx / det;
289 // compute alpha, a triangular square root
290 // of W (i.e. a Cholesky factor)
291 alpha(0, 0) = sqrt(W(0, 0));
292 // checked that alpha*alphaT = transW
293 alpha(1, 0) = W(0, 1) / alpha(0, 0);
294 alpha(1, 1) = 1. / sqrt(det * W(0, 0));
295 alpha(0, 1) = 0;
296 indices[0] = fs.getIndexInMatrix();
297 indices[1] = fs.getIndexInMatrix() + 1;
298 unsigned npar_tot = 2;
299 /* TODO: account here for proper motions in the reference
300 catalog. We can code the effect and set the value to 0. Most
301 (all?) catalogs do not even come with a reference epoch. Gaia
302 will change that. When refraction enters into the game, one should
303 pay attention to the orientation of the frame */
304
305 /* The residual should be Proj(fs)-Proj(*rs) in order to be consistent
306 with the measurement terms. Since P(fs) = 0, we have: */
307 res[0] = -rsProj.x;
308 res[1] = -rsProj.y;
309 halpha = H * alpha;
310 // grad = H*W*res
311 HW = H * W;
312 grad = HW * res;
313 // now feed in triplets and fullGrad
314 for (std::size_t ipar = 0; ipar < npar_tot; ++ipar) {
315 for (unsigned ic = 0; ic < 2; ++ic) {
316 double val = halpha(ipar, ic);
317 if (val == 0) continue;
318 tripletList.addTriplet(indices[ipar], kTriplets + ic, val);
319 }
320 fullGrad(indices[ipar]) += grad(ipar);
321 }
322 kTriplets += 2; // each measurement contributes 2 columns in the Jacobian
323 }
324 tripletList.setNextFreeIndex(kTriplets);
325}
326
327void AstrometryFit::accumulateStatImage(CcdImage const &ccdImage, Chi2Accumulator &accum) const {
328 /**********************************************************************/
331 /**********************************************************************/
332 /* Setup */
333 // 1 : get the Mapping's
334 const AstrometryMapping *mapping = _astrometryModel->getMapping(ccdImage);
335 // FittedStar is "observed" epoch, MeasuredStar is "baseline"
336 double deltaYears = _epoch - ccdImage.getEpoch();
337 // transformation from sky to TP
338 auto sky2TP = _astrometryModel->getSkyToTangentPlane(ccdImage);
339 // reserve matrix once for all measurements
340 Eigen::Matrix2Xd transW(2, 2);
341
342 auto &catalog = ccdImage.getCatalogForFit();
343 for (auto const &ms : catalog) {
344 if (!ms->isValid()) continue;
345 // tweak the measurement errors
346 FatPoint inPos = *ms;
347 tweakAstromMeasurementErrors(inPos, *ms, _posError);
348
349 FatPoint outPos;
350 // should *not* fill H if whatToFit excludes mapping parameters.
351 mapping->transformPosAndErrors(inPos, outPos);
352 double det = outPos.vx * outPos.vy - std::pow(outPos.vxy, 2);
353 if (det <= 0 || outPos.vx <= 0 || outPos.vy <= 0) {
354 LOGLS_WARN(_log, " Inconsistent measurement errors :drop measurement at "
355 << Point(*ms) << " in image " << ccdImage.getName());
356 continue;
357 }
358 transW(0, 0) = outPos.vy / det;
359 transW(1, 1) = outPos.vx / det;
360 transW(0, 1) = transW(1, 0) = -outPos.vxy / det;
361
362 std::shared_ptr<FittedStar const> const fs = ms->getFittedStar();
363 Point fittedStarInTP = transformFittedStar(*fs, *sky2TP, deltaYears);
364
365 Eigen::Vector2d res(fittedStarInTP.x - outPos.x, fittedStarInTP.y - outPos.y);
366 double chi2Val = res.transpose() * transW * res;
367
368 accum.addEntry(chi2Val, 2, ms);
369 } // end of loop on measurements
370}
371
373 for (auto const &ccdImage : ccdImageList) {
374 accumulateStatImage(*ccdImage, accum);
375 }
376}
377
379 /**********************************************************************/
382 /**********************************************************************/
383
384 /* If you wonder why we project here, read comments in
385 AstrometryFit::leastSquareDerivativesReference(TripletList &TList, Eigen::VectorXd &Rhs) */
386 FittedStarList &fittedStarList = _associations->fittedStarList;
388 for (auto const &fs : fittedStarList) {
389 auto rs = fs->getRefStar();
390 if (rs == nullptr) continue;
391 proj.setTangentPoint(*fs);
392 // fs projects to (0,0), no need to compute its transform.
393 FatPoint rsProj;
394 proj.transformPosAndErrors(*rs, rsProj);
395 // TO DO : account for proper motions.
396 double chi2 = computeProjectedRefStarChi2(rsProj);
397 accum.addEntry(chi2, 2, fs);
398 }
399}
400
402
404void AstrometryFit::getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const {
405 if (_fittingDistortions) {
406 const AstrometryMapping *mapping = _astrometryModel->getMapping(measuredStar.getCcdImage());
407 mapping->getMappingIndices(indices);
408 }
409 std::shared_ptr<FittedStar const> const fs = measuredStar.getFittedStar();
410 Eigen::Index fsIndex = fs->getIndexInMatrix();
411 if (_fittingPos) {
412 indices.push_back(fsIndex);
413 indices.push_back(fsIndex + 1);
414 }
415 // For securing the outlier removal, the next block is just useless
416 if (_fittingPM) {
417 for (std::size_t k = 0; k < 2; ++k) indices.push_back(fsIndex + 2 + k);
418 }
419 /* Should not put the index of refaction stuff or we will not be
420 able to remove more than 1 star at a time. */
421}
422
424 _whatToFit = whatToFit;
425 LOGLS_INFO(_log, "assignIndices: Now fitting " << whatToFit);
426 _fittingDistortions = (_whatToFit.find("Distortions") != std::string::npos);
427 _fittingPos = (_whatToFit.find("Positions") != std::string::npos);
428 _fittingPM = (_whatToFit.find("PM") != std::string::npos);
429 // When entering here, we assume that whatToFit has already been interpreted.
430
431 _nModelParams = 0;
432 if (_fittingDistortions) _nModelParams = _astrometryModel->assignIndices(_whatToFit, 0);
434
435 if (_fittingPos) {
436 FittedStarList &fittedStarList = _associations->fittedStarList;
437 for (auto &fittedStar : fittedStarList) {
438 // the parameter layout here is used also
439 // - when filling the derivatives
440 // - when updating (offsetParams())
441 // - in GetMeasuredStarIndices
442 fittedStar->setIndexInMatrix(ipar);
443 ipar += 2;
444 // TODO: left as reference for when we implement PM fitting
445 // if ((_fittingPM)&fittedStar->mightMove) ipar += 2;
446 }
447 }
449 _nTotal = ipar;
450 LOGLS_DEBUG(_log, "nParameters total: " << _nTotal << " model: " << _nModelParams
451 << " values: " << _nStarParams);
452}
453
454void AstrometryFit::offsetParams(Eigen::VectorXd const &delta) {
455 if (delta.size() != _nTotal)
457 "AstrometryFit::offsetParams : the provided vector length is not compatible with "
458 "the current whatToFit setting");
459 if (_fittingDistortions) _astrometryModel->offsetParams(delta);
460
461 if (_fittingPos) {
462 FittedStarList &fittedStarList = _associations->fittedStarList;
463 for (auto const &i : fittedStarList) {
464 FittedStar &fs = *i;
465 // the parameter layout here is used also
466 // - when filling the derivatives
467 // - when assigning indices (assignIndices())
468 Eigen::Index index = fs.getIndexInMatrix();
469 fs.x += delta(index);
470 fs.y += delta(index + 1);
471 // TODO: left as reference for when we implement PM fitting
472 // if ((_fittingPM)&fs.mightMove) {
473 // fs.pmx += delta(index + 2);
474 // fs.pmy += delta(index + 3);
475 // }
476 }
477 }
478}
479
481 const char *what2fit[] = {"Positions", "Distortions", "Positions Distortions"};
482 // DEBUG
483 for (unsigned k = 0; k < sizeof(what2fit) / sizeof(what2fit[0]); ++k) {
484 assignIndices(what2fit[k]);
485 TripletList tripletList(10000);
486 Eigen::VectorXd grad(_nTotal);
487 grad.setZero();
488 leastSquareDerivatives(tripletList, grad);
489 SparseMatrixD jacobian(_nTotal, tripletList.getNextFreeIndex());
490 jacobian.setFromTriplets(tripletList.begin(), tripletList.end());
491 SparseMatrixD hessian = jacobian * jacobian.transpose();
492 LOGLS_DEBUG(_log, "npar : " << _nTotal << ' ' << _nModelParams);
493 }
494}
495
497 std::ofstream ofile(filename.c_str());
498 std::string separator = "\t";
499
500 ofile << "id" << separator << "xccd" << separator << "yccd " << separator;
501 ofile << "rx" << separator << "ry" << separator;
502 ofile << "xtp" << separator << "ytp" << separator;
503 ofile << "mag" << separator << "deltaYears" << separator;
504 ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
505 ofile << "xtpi" << separator << "ytpi" << separator;
506 ofile << "rxi" << separator << "ryi" << separator;
507 ofile << "fsindex" << separator;
508 ofile << "ra" << separator << "dec" << separator;
509 ofile << "chi2" << separator << "nm" << separator;
510 ofile << "detector" << separator << "visit" << std::endl;
511
512 ofile << "id in source catalog" << separator << "coordinates in CCD (pixels)" << separator << separator;
513 ofile << "residual on TP (degrees)" << separator << separator;
514 ofile << "transformed coordinate in TP (degrees)" << separator << separator;
515 ofile << "rough magnitude" << separator << "Julian epoch year delta from fit epoch" << separator;
516 ofile << "transformed measurement uncertainty (degrees)" << separator << separator << separator;
517 ofile << "as-read position on TP (degrees)" << separator << separator;
518 ofile << "as-read residual on TP (degrees)" << separator << separator;
519 ofile << "unique index of the fittedStar" << separator;
520 ofile << "on sky position of fittedStar" << separator << separator;
521 ofile << "contribution to chi2 (2D dofs)" << separator << "number of measurements of this fittedStar"
522 << separator;
523 ofile << "detector id" << separator << "visit id" << std::endl;
524
525 const CcdImageList &ccdImageList = _associations->getCcdImageList();
526 for (auto const &ccdImage : ccdImageList) {
527 const MeasuredStarList &cat = ccdImage->getCatalogForFit();
528 const AstrometryMapping *mapping = _astrometryModel->getMapping(*ccdImage);
529 const auto readTransform = ccdImage->getReadWcs();
530 // FittedStar is "observed" epoch, MeasuredStar is "baseline"
531 double deltaYears = _epoch - ccdImage->getEpoch();
532 for (auto const &ms : cat) {
533 if (!ms->isValid()) continue;
534 FatPoint tpPos;
535 FatPoint inPos = *ms;
536 tweakAstromMeasurementErrors(inPos, *ms, _posError);
537 mapping->transformPosAndErrors(inPos, tpPos);
538 auto sky2TP = _astrometryModel->getSkyToTangentPlane(*ccdImage);
539 const std::unique_ptr<AstrometryTransform> readPixToTangentPlane =
540 compose(*sky2TP, *readTransform);
541 FatPoint inputTpPos = readPixToTangentPlane->apply(inPos);
542 std::shared_ptr<FittedStar const> const fs = ms->getFittedStar();
543
544 Point fittedStarInTP = transformFittedStar(*fs, *sky2TP, deltaYears);
545 Point res = tpPos - fittedStarInTP;
546 Point inputRes = inputTpPos - fittedStarInTP;
547 double det = tpPos.vx * tpPos.vy - std::pow(tpPos.vxy, 2);
548 double wxx = tpPos.vy / det;
549 double wyy = tpPos.vx / det;
550 double wxy = -tpPos.vxy / det;
551 double chi2 = wxx * res.x * res.x + wyy * res.y * res.y + 2 * wxy * res.x * res.y;
552
553 ofile << std::setprecision(9);
554 ofile << ms->getId() << separator << ms->x << separator << ms->y << separator;
555 ofile << res.x << separator << res.y << separator;
556 ofile << tpPos.x << separator << tpPos.y << separator;
557 ofile << fs->getMag() << separator << deltaYears << separator;
558 ofile << tpPos.vx << separator << tpPos.vy << separator << tpPos.vxy << separator;
559 ofile << inputTpPos.x << separator << inputTpPos.y << separator;
560 ofile << inputRes.x << separator << inputRes.y << separator;
561 ofile << fs->getIndexInMatrix() << separator;
562 ofile << fs->x << separator << fs->y << separator;
563 ofile << chi2 << separator << fs->getMeasurementCount() << separator;
564 ofile << ccdImage->getCcdId() << separator << ccdImage->getVisit() << std::endl;
565 } // loop on measurements in image
566 } // loop on images
567}
568
570 std::ofstream ofile(filename.c_str());
571 std::string separator = "\t";
572
573 ofile << "ra" << separator << "dec " << separator;
574 ofile << "rx" << separator << "ry" << separator;
575 ofile << "mag" << separator;
576 ofile << "xErr" << separator << "yErr" << separator << "xyCov" << separator;
577 ofile << "fsindex" << separator;
578 ofile << "chi2" << separator << "nm" << std::endl;
579
580 ofile << "coordinates of fittedStar (degrees)" << separator << separator;
581 ofile << "residual on TP (degrees)" << separator << separator;
582 ofile << "magnitude" << separator;
583 ofile << "refStar transformed measurement uncertainty (degrees)" << separator << separator << separator;
584 ofile << "unique index of the fittedStar" << separator;
585 ofile << "refStar contribution to chi2 (2D dofs)" << separator
586 << "number of measurements of this FittedStar" << std::endl;
587
588 // The following loop is heavily inspired from AstrometryFit::computeChi2()
589 const FittedStarList &fittedStarList = _associations->fittedStarList;
591 for (auto const &i : fittedStarList) {
592 const FittedStar &fs = *i;
593 auto rs = fs.getRefStar();
594 if (rs == nullptr) continue;
595 proj.setTangentPoint(fs);
596 // fs projects to (0,0), no need to compute its transform.
597 FatPoint rsProj;
598 proj.transformPosAndErrors(*rs, rsProj);
599 double chi2 = computeProjectedRefStarChi2(rsProj);
600
601 ofile << std::setprecision(9);
602 ofile << fs.x << separator << fs.y << separator;
603 ofile << rsProj.x << separator << rsProj.y << separator;
604 ofile << fs.getMag() << separator;
605 ofile << rsProj.vx << separator << rsProj.vy << separator << rsProj.vxy << separator;
606 ofile << fs.getIndexInMatrix() << separator;
607 ofile << chi2 << separator << fs.getMeasurementCount() << std::endl;
608 } // loop on FittedStars
609}
610} // namespace jointcal
611} // namespace lsst
Eigen::SparseMatrix< double, 0, Eigen::Index > SparseMatrixD
Definition Eigenstuff.h:35
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition Exception.h:48
LSST DM logging module built on log4cxx.
#define LOGLS_WARN(logger, message)
Log a warn-level message using an iostream-based interface.
Definition Log.h:659
#define LOG_GET(logger)
Returns a Log object associated with logger.
Definition Log.h:75
#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 at(T... args)
T begin(T... args)
T c_str(T... args)
void checkStuff()
DEBUGGING routine.
void leastSquareDerivativesMeasurement(CcdImage const &ccdImage, TripletList &tripletList, Eigen::VectorXd &grad, MeasuredStarList const *msList=nullptr) const override
Compute the derivatives of the measured stars and model for one CcdImage.
void accumulateStatImageList(CcdImageList const &ccdImageList, Chi2Accumulator &accum) const override
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for measurements.
void accumulateStatRefStars(Chi2Accumulator &accum) const override
Compute the chi2 (per star or total, depending on which Chi2Accumulator is used) for RefStars.
void saveChi2RefContributions(std::string const &filename) const override
Save a CSV file containing residuals of reference terms.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the requested quantities.
AstrometryFit(std::shared_ptr< Associations > associations, std::shared_ptr< AstrometryModel > astrometryModel, double posError)
this is the only constructor
void assignIndices(std::string const &whatToFit) override
Set parameters to fit and assign indices in the big matrix.
void leastSquareDerivativesReference(FittedStarList const &fittedStarList, TripletList &tripletList, Eigen::VectorXd &grad) const override
Compute the derivatives of the reference terms.
void getIndicesOfMeasuredStar(MeasuredStar const &measuredStar, IndexVector &indices) const override
this routine is to be used only in the framework of outlier removal
void saveChi2MeasContributions(std::string const &filename) const override
Save a CSV file containing residuals of measurement terms.
virtual class needed in the abstraction of the distortion model
virtual void transformPosAndErrors(FatPoint const &where, FatPoint &outPoint) const =0
The same as above but without the parameter derivatives (used to evaluate chi^2)
virtual void computeTransformAndDerivatives(FatPoint const &where, FatPoint &outPoint, Eigen::MatrixX2d &H) const =0
Actually applies the AstrometryMapping and evaluates the derivatives w.r.t the fitted parameters.
virtual void getMappingIndices(IndexVector &indices) const =0
Sets how this set of parameters (of length Npar()) map into the "grand" fit Expects that indices has ...
virtual std::size_t getNpar() const =0
Number of parameters in total.
a virtual (interface) class for geometric transformations.
virtual void apply(double xIn, double yIn, double &xOut, double &yOut) const =0
virtual void computeDerivative(Point const &where, AstrometryTransformLinear &derivative, double step=0.01) const
Computes the local Derivative of a transform, w.r.t.
implements the linear transformations (6 real coefficients).
double getMag() const
Definition BaseStar.h:105
Handler of an actual image from a single CCD.
Definition CcdImage.h:64
CcdIdType getCcdId() const
returns ccd ID
Definition CcdImage.h:145
VisitIdType getVisit() const
returns visit ID
Definition CcdImage.h:148
double getEpoch() const
Definition CcdImage.h:157
std::shared_ptr< AstrometryTransform > const getReadWcs() const
the wcs read in the header. NOT updated when fitting.
Definition CcdImage.h:187
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition CcdImage.h:94
std::string getName() const
Return the _name that identifies this ccdImage.
Definition CcdImage.h:79
Base class for Chi2Statistic and Chi2List, to allow addEntry inside Fitter for either class.
Definition Chi2.h:44
virtual void addEntry(double inc, std::size_t dof, std::shared_ptr< BaseStar > star)=0
A Point with uncertainties.
Definition FatPoint.h:34
FittedStars are objects whose position or flux is going to be fitted, and which come from the associa...
Definition FittedStar.h:50
void setIndexInMatrix(Eigen::Index const index)
index is a value that a fit can set and reread....
Definition FittedStar.h:96
const RefStar * getRefStar() const
Get the astrometric reference star associated with this star.
Definition FittedStar.h:105
int getMeasurementCount() const
The number of MeasuredStars currently associated with this FittedStar.
Definition FittedStar.h:88
Eigen::Index getIndexInMatrix() const
Definition FittedStar.h:99
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition FittedStar.h:116
Base class for fitters.
Definition FitterBase.h:55
void leastSquareDerivatives(TripletList &tripletList, Eigen::VectorXd &grad) const
Evaluates the chI^2 derivatives (Jacobian and gradient) for the current whatToFit setting.
std::shared_ptr< Associations > _associations
Definition FitterBase.h:165
Sources measured on images.
bool isValid() const
Fits may use that to discard outliers.
CcdImage const & getCcdImage() const
std::shared_ptr< FittedStar > getFittedStar() const
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
A point in a plane.
Definition Point.h:37
double x
coordinate
Definition Point.h:42
Point applyProperMotion(const Point &star, double timeDeltaYears) const
Apply proper motion correction to the input star, returning a star with PM-corrected coordinates and ...
Definition RefStar.cc:34
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane)
void transformPosAndErrors(const FatPoint &in, FatPoint &out) const
transform with analytical derivatives
void setTangentPoint(Point const &tangentPoint)
Resets the projection (or tangent) point.
Eigen::Index getNextFreeIndex() const
Definition Tripletlist.h:47
void setNextFreeIndex(Eigen::Index index)
Definition Tripletlist.h:49
void addTriplet(Eigen::Index i, Eigen::Index j, double val)
Definition Tripletlist.h:43
Reports invalid arguments.
Definition Runtime.h:66
T end(T... args)
T endl(T... args)
T find(T... args)
T front(T... args)
T move(T... args)
std::unique_ptr< AstrometryTransform > compose(AstrometryTransform const &left, AstrometryTransform const &right)
Returns a pointer to a composition of transforms, representing left(right()).
STL namespace.
T pow(T... args)
T push_back(T... args)
T setprecision(T... args)
T sqrt(T... args)
ImageT val
Definition CR.cc:146