lsst.jointcal  16.0-24-gd8faad8+3
Associations.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 <cmath>
26 #include <iostream>
27 #include <limits>
28 #include <sstream>
29 
30 #include "lsst/log/Log.h"
32 #include "lsst/jointcal/CcdImage.h"
35 #include "lsst/jointcal/Frame.h"
36 #include "lsst/jointcal/FatPoint.h"
37 #include "lsst/jointcal/Gtransfo.h"
39 
40 #include "lsst/afw/image/Image.h"
43 
44 #include "lsst/pex/exceptions.h"
45 #include "lsst/afw/geom/Box.h"
46 #include "lsst/afw/geom/Point.h"
47 #include "lsst/afw/image/Calib.h"
48 
49 namespace jointcal = lsst::jointcal;
50 
51 namespace {
52 LOG_LOGGER _log = LOG_GET("jointcal.Associations");
53 }
54 
55 // TODO: Remove this once RFC-356 is implemented and all refcats give fluxes in Maggies.
56 const double JanskyToMaggy = 3631.0;
57 
58 namespace lsst {
59 namespace jointcal {
60 
64  lsst::afw::geom::Box2I const &bbox, std::string const &filter,
66  std::shared_ptr<afw::cameraGeom::Detector> detector, int visit, int ccd,
67  lsst::jointcal::JointcalControl const &control) {
68  auto ccdImage = std::make_shared<CcdImage>(catalog, wcs, visitInfo, bbox, filter, photoCalib, detector,
69  visit, ccd, control.sourceFluxField);
70  ccdImageList.push_back(ccdImage);
71  LOGLS_DEBUG(_log, "Catalog " << ccdImage->getName() << " has " << ccdImage->getWholeCatalog().size()
72  << " objects.");
73 }
74 
77  centers.reserve(ccdImageList.size());
78  for (auto const &ccdImage : ccdImageList) {
79  centers.push_back(ccdImage->getBoresightRaDec());
80  }
81  auto commonTangentPoint = afw::geom::averageSpherePoint(centers);
82  LOGLS_DEBUG(_log, "Using common tangent point: " << commonTangentPoint.getPosition(afw::geom::degrees));
83  setCommonTangentPoint(commonTangentPoint.getPosition(afw::geom::degrees));
84 }
85 
86 void Associations::setCommonTangentPoint(lsst::afw::geom::Point2D const &commonTangentPoint) {
87  _commonTangentPoint = Point(commonTangentPoint.getX(), commonTangentPoint.getY()); // a jointcal::Point
88  for (auto &ccdImage : ccdImageList) ccdImage->setCommonTangentPoint(_commonTangentPoint);
89 }
90 
91 void Associations::associateCatalogs(const double matchCutInArcSec, const bool useFittedList,
92  const bool enlargeFittedList) {
93  // clear reference stars
94  refStarList.clear();
95 
96  // clear measurement counts and associations to refstars, but keep fittedStars themselves.
97  for (auto &item : fittedStarList) {
98  item->clearBeforeAssoc();
99  }
100  // clear fitted stars
101  if (!useFittedList) fittedStarList.clear();
102 
103  for (auto &ccdImage : ccdImageList) {
104  const Gtransfo *toCommonTangentPlane = ccdImage->getPix2CommonTangentPlane();
105 
106  // Clear the catalog to fit and copy the whole catalog into it.
107  // This allows reassociating from scratch after a fit.
108  ccdImage->resetCatalogForFit();
109  MeasuredStarList &catalog = ccdImage->getCatalogForFit();
110 
111  // Associate with previous lists.
112  /* To speed up the match (more precisely the contruction of the FastFinder), select in the
113  fittedStarList the objects that are within reach of the current ccdImage */
114  Frame ccdImageFrameCPT = toCommonTangentPlane->apply(ccdImage->getImageFrame(), false);
115  ccdImageFrameCPT = ccdImageFrameCPT.rescale(1.10); // add 10 % margin.
116  // We cannot use FittedStarList::ExtractInFrame, because it does an actual copy, which we don't want
117  // here: we want the pointers in the StarMatch to refer to fittedStarList elements.
118  FittedStarList toMatch;
119 
120  for (auto const &fittedStar : fittedStarList) {
121  if (ccdImageFrameCPT.inFrame(*fittedStar)) {
122  toMatch.push_back(fittedStar);
123  }
124  }
125 
126  // divide by 3600 because coordinates in CTP are in degrees.
127  auto starMatchList = listMatchCollect(Measured2Base(catalog), Fitted2Base(toMatch),
128  toCommonTangentPlane, matchCutInArcSec / 3600.);
129 
130  /* should check what this removeAmbiguities does... */
131  LOGLS_DEBUG(_log, "Measured-to-Fitted matches before removing ambiguities " << starMatchList->size());
132  starMatchList->removeAmbiguities(*toCommonTangentPlane);
133  LOGLS_DEBUG(_log, "Measured-to-Fitted matches after removing ambiguities " << starMatchList->size());
134 
135  // Associate MeasuredStar -> FittedStar using the surviving matches.
136 
137  int matchedCount = 0;
138  for (auto const &starMatch : *starMatchList) {
139  auto bs = starMatch.s1;
140  auto ms_const = std::dynamic_pointer_cast<const MeasuredStar>(bs);
141  auto ms = std::const_pointer_cast<MeasuredStar>(ms_const);
142  auto bs2 = starMatch.s2;
143  auto fs_const = std::dynamic_pointer_cast<const FittedStar>(bs2);
144  auto fs = std::const_pointer_cast<FittedStar>(fs_const);
145  ms->setFittedStar(fs);
146  matchedCount++;
147  }
148  LOGLS_INFO(_log, "Matched " << matchedCount << " objects in " << ccdImage->getName());
149 
150  // add unmatched objets to FittedStarList
151  int unMatchedCount = 0;
152  for (auto const &mstar : catalog) {
153  // to check if it was matched, just check if it has a fittedStar Pointer assigned
154  if (mstar->getFittedStar()) continue;
155  if (enlargeFittedList) {
156  auto fs = std::make_shared<FittedStar>(*mstar);
157  // transform coordinates to CommonTangentPlane
158  toCommonTangentPlane->transformPosAndErrors(*fs, *fs);
159  fittedStarList.push_back(fs);
160  mstar->setFittedStar(fs);
161  }
162  unMatchedCount++;
163  }
164  LOGLS_INFO(_log, "Unmatched objects: " << unMatchedCount);
165  } // end of loop on CcdImages
166 
167  // !!!!!!!!!!!!!!!!!
168  // TODO: DO WE REALLY NEED THIS???
169  // Why do we need to do this, instead of directly computing them in normalizeFittedStars?
170  // What makes the magnitudes special here?
171  // !!!!!!!!!!!!!!!!!
172  // assignMags();
173 }
174 
175 void Associations::collectRefStars(afw::table::SimpleCatalog &refCat, afw::geom::Angle matchCut,
176  std::string const &fluxField,
177  std::map<std::string, std::vector<double>> const &refFluxMap,
178  std::map<std::string, std::vector<double>> const &refFluxErrMap,
179  bool rejectBadFluxes) {
180  if (refCat.size() == 0) {
182  " reference catalog is empty : stop here "));
183  }
184 
185  afw::table::CoordKey coordKey = refCat.getSchema()["coord"];
186  auto fluxKey = refCat.getSchema().find<double>(fluxField).key;
187  // Don't blow up if the reference catalog doesn't contain errors.
188  afw::table::Key<double> fluxErrKey;
189  try {
190  fluxErrKey = refCat.getSchema().find<double>(fluxField + "Err").key;
191  } catch (pex::exceptions::NotFoundError &) {
192  LOGLS_WARN(_log, "Flux error field ("
193  << fluxField << "Err"
194  << ") not found in reference catalog. Not using ref flux errors.");
195  }
196  _filterMap.clear();
197  _filterMap.reserve(refFluxMap.size());
198  size_t nFilters = 0;
199  for (auto const &filter : refFluxMap) {
200  _filterMap[filter.first] = nFilters;
201  nFilters++;
202  }
203 
204  refStarList.clear();
205  for (size_t i = 0; i < refCat.size(); i++) {
206  auto const &record = refCat.get(i);
207 
208  auto coord = record->get(coordKey);
209  double defaultFlux = record->get(fluxKey) / JanskyToMaggy;
210  double defaultFluxErr;
211  if (fluxErrKey.isValid()) {
212  defaultFluxErr = record->get(fluxErrKey) / JanskyToMaggy;
213  } else {
214  defaultFluxErr = std::numeric_limits<double>::quiet_NaN();
215  }
216  std::vector<double> fluxList(nFilters);
217  std::vector<double> fluxErrList(nFilters);
218  for (auto const &filter : _filterMap) {
219  fluxList[filter.second] = refFluxMap.at(filter.first).at(i) / JanskyToMaggy;
220  fluxErrList[filter.second] = refFluxErrMap.at(filter.first).at(i) / JanskyToMaggy;
221  }
222  double ra = lsst::afw::geom::radToDeg(coord.getLongitude());
223  double dec = lsst::afw::geom::radToDeg(coord.getLatitude());
224  auto star = std::make_shared<RefStar>(ra, dec, defaultFlux, defaultFluxErr, fluxList, fluxErrList);
225 
226  // TODO DM-10826: RefCats aren't guaranteed to have position errors.
227  // TODO: Need to devise a way to check whether the refCat has position errors
228  // TODO: and use them instead, if available.
229  // cook up errors: 100 mas per cooordinate
230  star->vx = std::pow(0.1 / 3600 / cos(coord.getLatitude()), 2);
231  star->vy = std::pow(0.1 / 3600, 2);
232  star->vxy = 0.;
233 
234  // Reject sources with non-finite fluxes and flux errors, and fluxErr=0 (which gives chi2=inf).
235  if (rejectBadFluxes &&
236  (!std::isfinite(defaultFlux) || !std::isfinite(defaultFluxErr) || defaultFluxErr == 0))
237  continue;
238  refStarList.push_back(star);
239  }
240 
241  // project on CTP (i.e. RaDec2CTP), in degrees
242  GtransfoLin identity;
243  TanRaDec2Pix raDec2CTP(identity, _commonTangentPoint);
244 
245  associateRefStars(matchCut.asArcseconds(), &raDec2CTP);
246 }
247 
248 const lsst::afw::geom::Box2D Associations::getRaDecBBox() {
249  // compute the frame on the CTP that contains all input images
250  Frame tangentPlaneFrame;
251 
252  for (auto const &ccdImage : ccdImageList) {
253  Frame CTPFrame = ccdImage->getPix2CommonTangentPlane()->apply(ccdImage->getImageFrame(), false);
254  if (tangentPlaneFrame.getArea() == 0)
255  tangentPlaneFrame = CTPFrame;
256  else
257  tangentPlaneFrame += CTPFrame;
258  }
259 
260  // convert tangent plane coordinates to RaDec:
261  GtransfoLin identity;
262  TanPix2RaDec CTP2RaDec(identity, _commonTangentPoint);
263  Frame raDecFrame = CTP2RaDec.apply(tangentPlaneFrame, false);
264 
265  lsst::afw::geom::Point<double> min(raDecFrame.xMin, raDecFrame.yMin);
266  lsst::afw::geom::Point<double> max(raDecFrame.xMax, raDecFrame.yMax);
267  lsst::afw::geom::Box2D box(min, max);
268 
269  return box;
270 }
271 
272 void Associations::associateRefStars(double matchCutInArcSec, const Gtransfo *gtransfo) {
273  // associate with FittedStars
274  // 3600 because coordinates are in degrees (in CTP).
275  auto starMatchList = listMatchCollect(Ref2Base(refStarList), Fitted2Base(fittedStarList), gtransfo,
276  matchCutInArcSec / 3600.);
277 
278  LOGLS_DEBUG(_log, "Refcat matches before removing ambiguities " << starMatchList->size());
279  starMatchList->removeAmbiguities(*gtransfo);
280  LOGLS_DEBUG(_log, "Refcat matches after removing ambiguities " << starMatchList->size());
281 
282  // actually associate things
283  for (auto const &starMatch : *starMatchList) {
284  const BaseStar &bs = *starMatch.s1;
285  const RefStar &rs_const = dynamic_cast<const RefStar &>(bs);
286  RefStar &rs = const_cast<RefStar &>(rs_const);
287  const BaseStar &bs2 = *starMatch.s2;
288  const FittedStar &fs_const = dynamic_cast<const FittedStar &>(bs2);
289  FittedStar &fs = const_cast<FittedStar &>(fs_const);
290  // rs->setFittedStar(*fs);
291  fs.setRefStar(&rs);
292  }
293 
294  LOGLS_INFO(_log,
295  "Associated " << starMatchList->size() << " reference stars among " << refStarList.size());
296 }
297 
298 void Associations::prepareFittedStars(int minMeasurements) {
299  selectFittedStars(minMeasurements);
300  normalizeFittedStars();
301 }
302 
303 void Associations::selectFittedStars(int minMeasurements) {
304  LOGLS_INFO(_log, "Fitted stars before measurement # cut: " << fittedStarList.size());
305 
306  // first pass: remove objects that have less than a certain number of measurements.
307  for (auto const &ccdImage : ccdImageList) {
308  MeasuredStarList &catalog = ccdImage->getCatalogForFit();
309  // Iteration happens internal to the loop, as we may delete measuredStars from catalog.
310  for (MeasuredStarIterator mi = catalog.begin(); mi != catalog.end();) {
311  MeasuredStar &mstar = **mi;
312 
313  auto fittedStar = mstar.getFittedStar();
314  // measuredStar has no fittedStar: move on.
315  if (fittedStar == nullptr) {
316  ++mi;
317  continue;
318  }
319 
320  // keep FittedStars which either have a minimum number of
321  // measurements, or are matched to a RefStar
322  if (!fittedStar->getRefStar() && fittedStar->getMeasurementCount() < minMeasurements) {
323  fittedStar->getMeasurementCount()--;
324  mi = catalog.erase(mi); // mi now points to the next measuredStar.
325  } else {
326  ++mi;
327  }
328  } // end loop on objects in catalog
329  } // end loop on catalogs
330 
331  // now FittedStars with less than minMeasurements should have zero measurementCount.
332  for (FittedStarIterator fi = fittedStarList.begin(); fi != fittedStarList.end();) {
333  if ((*fi)->getMeasurementCount() == 0) {
334  fi = fittedStarList.erase(fi);
335  } else {
336  ++fi;
337  }
338  }
339 
340  LOGLS_INFO(_log, "Fitted stars after measurement # cut: " << fittedStarList.size());
341 }
342 
343 void Associations::normalizeFittedStars() const {
344  // Clear positions in order to take the average of the measuredStars.
345  for (auto &fittedStar : fittedStarList) {
346  fittedStar->x = 0.0;
347  fittedStar->y = 0.0;
348  fittedStar->setFlux(0.0);
349  fittedStar->getMag() = 0.0;
350  }
351 
352  // Iterate over measuredStars to add their values into their fittedStars
353  for (auto const &ccdImage : ccdImageList) {
354  const Gtransfo *toCommonTangentPlane = ccdImage->getPix2CommonTangentPlane();
355  MeasuredStarList &catalog = ccdImage->getCatalogForFit();
356  for (auto &mi : catalog) {
357  auto fittedStar = mi->getFittedStar();
358  if (fittedStar == nullptr)
359  throw(LSST_EXCEPT(
361  "All measuredStars must have a fittedStar: did you call selectFittedStars()?"));
362  auto point = toCommonTangentPlane->apply(*mi);
363  fittedStar->x += point.x;
364  fittedStar->y += point.y;
365  fittedStar->getFlux() += mi->getFlux();
366  }
367  }
368 
369  for (auto &fi : fittedStarList) {
370  auto measurementCount = fi->getMeasurementCount();
371  fi->x /= measurementCount;
372  fi->y /= measurementCount;
373  fi->getFlux() /= measurementCount;
374  fi->getMag() = magFromFlux(fi->getFlux());
375  }
376 }
377 
378 void Associations::assignMags() {
379  for (auto const &ccdImage : ccdImageList) {
380  MeasuredStarList &catalog = ccdImage->getCatalogForFit();
381  for (auto const &mstar : catalog) {
382  auto fstar = mstar->getFittedStar();
383  if (!fstar) continue;
384  fstar->addMagMeasurement(mstar->getMag(), mstar->getMagWeight());
385  }
386  }
387 }
388 
390  // By default, Associations::fittedStarList is expressed on the Associations::commonTangentPlane.
391  // For AstrometryFit, we need it on the sky.
392  if (!fittedStarList.inTangentPlaneCoordinates) {
393  LOGLS_WARN(_log,
394  "DeprojectFittedStars: Fitted stars are already in sidereal coordinates, nothing done ");
395  return;
396  }
397 
398  TanPix2RaDec ctp2Sky(GtransfoLin(), getCommonTangentPoint());
399  fittedStarList.applyTransfo(ctp2Sky);
400  fittedStarList.inTangentPlaneCoordinates = false;
401 }
402 
404  return std::count_if(ccdImageList.begin(), ccdImageList.end(), [](std::shared_ptr<CcdImage> const &item) {
405  return item->getCatalogForFit().size() > 0;
406  });
407 }
408 
410  size_t count = 0;
411  for (auto const &fittedStar : fittedStarList) {
412  if ((fittedStar != nullptr) & (fittedStar->getRefStar() != nullptr)) count++;
413  }
414  return count;
415 }
416 
417 #ifdef TODO
418 void Associations::collectMCStars(int realization) {
419  CcdImageIterator I;
420  StarMatchIterator smI;
421 
422  for (I = ccdImageList.begin(); I != ccdImageList.end(); I++) {
423  CcdImage &ccdImage = **I;
424  string dbimdir = ccdImage.Dir();
425  string mctruth = dbimdir + "/mc/mctruth.list";
426 
427  if (realization >= 0) {
428  stringstream sstrm;
429  sstrm << dbimdir << "/mc/mctruth_" << realization << ".list";
430  mctruth = sstrm.str();
431  }
432 
433  GtransfoIdentity gti;
434  MeasuredStarList &catalog = ccdImage.getCatalogForFit();
435 
436  // BaseStarWithErrorList mctruthlist(mctruth);
437  DicStarList mctruthlist(mctruth);
438  auto starMatchList =
439  listMatchCollect(Measured2Base(catalog), Dic2Base(mctruthlist), &gti, 1. /* pixel ? */);
440  if (starMatchList)
441  for (smI = starMatchList->begin(); smI != starMatchList->end(); smI++) {
442  StarMatch &sm = *smI;
443  BaseStar *bs = sm.s1;
444  MeasuredStar *mstar = dynamic_cast<MeasuredStar *>(bs);
445  bs = sm.s2;
446  DicStar *dstar = dynamic_cast<DicStar *>(bs);
447  std::unique_ptr<BaseStarWithError> mcstar(new BaseStarWithError(*bs));
448  mcstar->GetMCInfo().iflux = dstar->getval("iflux");
449  mcstar->GetMCInfo().tflux = dstar->getval("sflux");
450  /*
451  mstar->SetMCTruth(mcstar);
452  mstar->SetMCMeas(mcstar);
453  */
454  }
455  else
456  LOGLS_FATAL(_log, "CollectMCStars Unable to match MCTruth w/ catalog!");
457  }
458 }
459 
460 void Associations::setFittedStarColors(std::string dicStarListName, std::string color,
461  double matchCutArcSec) {
462  // decode color string in case it is x-y
463  size_t pos_minus = color.find('-');
464  bool compute_diff = (pos_minus != string::npos);
465  std::string c1, c2;
466  c1 = color.substr(0, pos_minus); // if pos_minus == npos, means "up to the end"
467  if (compute_diff) c2 = color.substr(pos_minus + 1, string::npos);
468  DicStarList cList(dicStarListName);
469  if (!cList.HasKey(c1))
470  throw(GastroException("Associations::SetFittedstarColors : " + dicStarListName +
471  " misses a key named \"" + c1 + "\""));
472  if (compute_diff && !cList.HasKey(c2))
473  throw(GastroException("Associations::SetFittedstarColors : " + dicStarListName +
474  " misses a key named \"" + c2 + "\""));
475  // we associate in some tangent plane. The reference catalog is expressed on the sky,
476  // but FittedStar's may be still in this tangent plane.
477  BaseStarList &l1 = (BaseStarList &)fittedStarList;
479  TanRaDec2Pix proj(GtransfoLin(), getCommonTangentPoint());
480  // project or not ?
481  Gtransfo *id_or_proj = &proj;
482  if (fittedStarList.inTangentPlaneCoordinates) id_or_proj = &id;
483  // The color List is to be projected:
484  TStarList projected_cList((BaseStarList &)cList, proj);
485  // Associate
486  auto starMatchList = listMatchCollect(Fitted2Base(fittedStarList), (const BaseStarList &)projected_cList,
487  id_or_proj, matchCutArcSec / 3600);
488 
489  LOGLS_INFO(_log, "Matched " << starMatchList->size() << '/' << fittedStarList.size()
490  << " FittedStars to color catalog");
491  // Evaluate and assign colors.
492  for (auto i = starMatchList->begin(); i != starMatchList->end(); ++i) {
493  BaseStar *s1 = i->s1;
494  FittedStar *fs = dynamic_cast<FittedStar *>(s1);
495  BaseStar *s2 = i->s2;
496  const TStar *ts = dynamic_cast<const TStar *>(s2);
497  const DicStar *ds = dynamic_cast<const DicStar *>(ts->get_original());
498  fs->color = ds->getval(c1);
499  if (compute_diff) fs->color -= ds->getval(c2);
500  }
501 }
502 
503 #endif /* TODO */
504 } // namespace jointcal
505 } // namespace lsst
#define LOGLS_WARN(logger, message)
Objects used as position anchors, typically USNO stars.
Definition: RefStar.h:39
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:414
double dec
int nCcdImagesValidForFit() const
return the number of CcdImages with non-empty catalogs to-be-fit.
MeasuredStarList::iterator MeasuredStarIterator
Definition: MeasuredStar.h:149
A hanger for star associations.
Definition: StarMatch.h:54
A point in a plane.
Definition: Point.h:36
double getArea() const
Definition: Frame.cc:118
the transformation that handles pix to sideral transfos (Gnomonic, possibly with polynomial distortio...
Definition: Gtransfo.h:600
void apply(const double xIn, const double yIn, double &xOut, double &yOut) const
Transform pixels to ICRS RA, Dec in degrees.
Definition: Gtransfo.cc:1351
double xMin
coordinate of boundary.
Definition: Frame.h:41
T end(T... args)
STL class.
void associateCatalogs(const double matchCutInArcsec=0, const bool useFittedList=false, const bool enlargeFittedList=true)
incrementaly builds a merged catalog of all image catalogs
Definition: Associations.cc:91
A list of MeasuredStar. They are usually filled in Associations::createCcdImage.
Definition: MeasuredStar.h:141
Frame rescale(const double factor) const
rescale it. The center does not move.
Definition: Frame.cc:110
STL class.
pairs of points
The base class for handling stars. Used by all matching routines.
Definition: BaseStar.h:50
T at(T... args)
void prepareFittedStars(int minMeasurements)
Set the color field of FittedStar &#39;s from a colored catalog.
size_t nFittedStarsWithAssociatedRefStar() const
Return the number of fittedStars that have an associated refStar.
T push_back(T... args)
void collectRefStars(afw::table::SimpleCatalog &refCat, afw::geom::Angle matchCut, std::string const &fluxField, RefFluxMapType const &refFluxMap=RefFluxMapType(), RefFluxMapType const &refFluxErrMap=RefFluxMapType(), bool rejectBadFluxes=false)
Collect stars from an external reference catalog and associate them with fittedStars.
std::lists of Stars.
Definition: StarList.h:58
rectangle with sides parallel to axes.
Definition: Frame.h:38
#define LOGLS_DEBUG(logger, message)
SchemaItem< T > find(std::string const &name) const
Class for a simple mapping implementing a generic Gtransfo.
void setRefStar(const RefStar *_refStar)
Set the astrometric reference star associated with this star.
Definition: FittedStar.cc:46
void computeCommonTangentPoint()
Sets a shared tangent point for all ccdImages, using the mean of the centers of all ccdImages...
Definition: Associations.cc:75
FittedStarList::iterator FittedStarIterator
Definition: FittedStar.h:131
int max
T erase(T... args)
A list of FittedStar s. Such a list is typically constructed by Associations.
Definition: FittedStar.h:122
Key< U > key
std::unique_ptr< StarMatchList > listMatchCollect(const BaseStarList &list1, const BaseStarList &list2, const Gtransfo *guess, const double maxDist)
assembles star matches.
Definition: ListMatch.cc:562
T str(T... args)
T isfinite(T... args)
virtual void transformPosAndErrors(const FatPoint &in, FatPoint &out) const
Definition: Gtransfo.cc:140
T cos(T... args)
objects measured on actual images.
Definition: MeasuredStar.h:42
T dynamic_pointer_cast(T... args)
::std::list< StarMatch >::iterator StarMatchIterator
Definition: StarMatch.h:133
void deprojectFittedStars()
Sends back the fitted stars coordinates on the sky FittedStarsList::inTangentPlaneCoordinates keeps t...
table::Key< int > id
T count_if(T... args)
#define LOGLS_INFO(logger, message)
void setCommonTangentPoint(lsst::afw::geom::Point2D const &commonTangentPoint)
Sets a shared tangent point for all ccdImages.
Definition: Associations.cc:86
bool inFrame(double x, double y) const
inside?
Definition: Frame.cc:120
T find(T... args)
A do-nothing transformation. It anyway has dummy routines to mimick a Gtransfo.
Definition: Gtransfo.h:219
#define LSST_EXCEPT(type,...)
void createCcdImage(afw::table::SourceCatalog &catalog, std::shared_ptr< lsst::afw::geom::SkyWcs > wcs, std::shared_ptr< lsst::afw::image::VisitInfo > visitInfo, lsst::afw::geom::Box2I const &bbox, std::string const &filter, std::shared_ptr< afw::image::PhotoCalib > photoCalib, std::shared_ptr< afw::cameraGeom::Detector > detector, int visit, int ccd, lsst::jointcal::JointcalControl const &control)
Create a ccdImage from an exposure catalog and metadata, and add it to the list.
Definition: Associations.cc:61
STL class.
STL class.
int min
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:676
T begin(T... args)
T pow(T... args)
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:94
Combinatorial searches for linear transformations to go from list1 to list2.
BaseStarList & Fitted2Base(FittedStarList &This)
Definition: FittedStar.cc:64
a virtual (interface) class for geometric transformations.
Definition: Gtransfo.h:65
std::shared_ptr< const BaseStar > s2
Definition: StarMatch.h:62
const lsst::afw::geom::Box2D getRaDecBBox()
BaseStarList & Measured2Base(MeasuredStarList &This)
Definition: MeasuredStar.cc:58
std::shared_ptr< RecordT > const get(size_type i) const
std::shared_ptr< const BaseStar > s1
Definition: StarMatch.h:62
T substr(T... args)
std::string sourceFluxField
"name of flux field in source catalog" ;
T quiet_NaN(T... args)
#define LOGLS_FATAL(logger, message)
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
#define LOG_GET(logger)
const double JanskyToMaggy
Definition: Associations.cc:56
std::shared_ptr< FittedStar > getFittedStar() const
Definition: MeasuredStar.h:108
The objects which have been measured several times.
Definition: FittedStar.h:60
virtual void apply(const double xIn, const double yIn, double &xOut, double &yOut) const =0
T reserve(T... args)
BaseStarList & Ref2Base(RefStarList &This)
Definition: RefStar.cc:34