lsst.jointcal  14.0-20-g39ceace
SimplePolyModel.cc
Go to the documentation of this file.
1 #include <string>
2 
3 #include "lsst/log/Log.h"
9 #include "lsst/pex/exceptions.h"
10 #include "lsst/jointcal/Gtransfo.h"
11 
12 // const int distortionDegree=3;
13 
14 namespace {
15 LOG_LOGGER _log = LOG_GET("jointcal.SimplePolyModel");
16 }
17 
18 namespace lsst {
19 namespace jointcal {
20 
21 // need a way to propagate the requested degree !
23  bool initFromWcs, unsigned nNotFit, unsigned degree)
24  : _sky2TP(projectionHandler)
25 
26 {
27  // from datacards (or default)
28  // unsigned degree = distortionDegree;
29  unsigned count = 0;
30 
31  for (auto i = ccdImageList.cbegin(); i != ccdImageList.cend(); ++i, ++count) {
32  const CcdImage &im = **i;
33  if (count < nNotFit) {
35  id->setIndex(-1); // non sense, because it has no parameters
36  _myMap[&im] = std::move(id);
37  } else
38  // Given how AssignIndices works, only the SimplePolyMapping's
39  // will actually be fitted, as nNotFit requests.
40  {
41  /* first check that there are enough measurements for the
42  requested polynomial degree */
43  size_t nObj = im.getCatalogForFit().size();
44  if (nObj == 0) {
45  LOGLS_WARN(_log, "Empty catalog from image: " << im.getName());
46  continue;
47  }
48  GtransfoPoly pol(degree);
49  if (pol.getDegree() > 0) // if not, it cannot be decreased
50  while (unsigned(pol.getNpar()) > 2 * nObj) pol.setDegree(pol.getDegree() - 1);
51  /* We have to center and normalize the coordinates so that
52  the fit matrix is not too ill-conditionned. Basically, x
53  and y in pixels are mapped to [-1,1]. When the
54  transformation of SimplePolyMapping transformation is
55  accessed, the combination of the normalization and the
56  fitted transformation is returned, so that the trick
57  remains hidden
58  */
59  const Frame &frame = im.getImageFrame();
60  GtransfoLin shiftAndNormalize = normalizeCoordinatesTransfo(frame);
61  if (initFromWcs) {
62  pol = GtransfoPoly(im.getPix2TangentPlane(), frame, degree);
63  pol = pol * shiftAndNormalize.invert();
64  }
65  _myMap[&im] =
66  std::unique_ptr<SimpleGtransfoMapping>(new SimplePolyMapping(shiftAndNormalize, pol));
67  }
68  }
69 }
70 
72  mapType::const_iterator i = _myMap.find(&ccdImage);
73  if (i == _myMap.cend())
75  "SimplePolyModel::GetMapping, never heard of CcdImage " + ccdImage.getName());
76  return (i->second.get());
77 }
78 
79 unsigned SimplePolyModel::assignIndices(unsigned firstIndex, std::string const &whatToFit) {
80  if (whatToFit.find("Distortions") == std::string::npos) {
81  LOGLS_ERROR(_log, "AssignIndices was called and Distortions is *not* in whatToFit.");
82  return 0;
83  }
84  unsigned index = firstIndex;
85  for (auto i = _myMap.begin(); i != _myMap.end(); ++i) {
86  SimplePolyMapping *p = dynamic_cast<SimplePolyMapping *>(&*(i->second));
87  if (!p) continue; // it should be GtransfoIdentity
88  p->setIndex(index);
89  index += p->getNpar();
90  }
91  return index;
92 }
93 
94 void SimplePolyModel::offsetParams(Eigen::VectorXd const &delta) {
95  for (auto &i : _myMap) {
96  auto mapping = i.second.get();
97  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
98  }
99 }
100 
102  for (auto i = _myMap.begin(); i != _myMap.end(); ++i) i->second->freezeErrorTransform();
103 }
104 
106  // return GetMapping(ccdImage)->Transfo(); // cannot do that
107  auto p = _myMap.find(&ccdImage);
108  if (p == _myMap.end())
110  "SimplePolyModel::getTransfo, never heard of CcdImage " + ccdImage.getName());
111  return p->second->getTransfo();
112 }
113 
115  const GtransfoPoly &pix2Tp = dynamic_cast<const GtransfoPoly &>(getTransfo(ccdImage));
116  const TanRaDec2Pix *proj = dynamic_cast<const TanRaDec2Pix *>(getSky2TP(ccdImage));
117  if (!proj) return nullptr;
118 
119  const GtransfoLin &projLinPart = proj->getLinPart(); // should be the identity, but who knows? So, let us
120  // incorporate it into the pix2TP part.
121  GtransfoPoly wcsPix2Tp = GtransfoPoly(projLinPart.invert()) * pix2Tp;
122 
123  // compute a decent approximation, if higher order corrections get ignored
124  GtransfoLin cdStuff = wcsPix2Tp.linearApproximation(ccdImage.getImageFrame().getCenter());
125 
126  // wcsPix2TP = cdStuff*sip , so
127  GtransfoPoly sip = GtransfoPoly(cdStuff.invert()) * wcsPix2Tp;
128  Point tangentPoint(proj->getTangentPoint());
129  return std::make_shared<TanSipPix2RaDec>(cdStuff, tangentPoint, &sip);
130 }
131 } // namespace jointcal
132 } // namespace lsst
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:292
virtual class needed in the abstraction of the distortion model
Definition: Mapping.h:15
A point in a plane.
Definition: Point.h:13
GtransfoLin invert() const
returns the inverse: T1 = T2.invert();
Definition: Gtransfo.cc:1089
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:50
Mapping implementation for a polynomial transformation.
const Gtransfo * getSky2TP(CcdImage const &ccdImage) const
the mapping of sky coordinates (i.e.
T cend(T... args)
Gtransfo const * getPix2TangentPlane() const
Definition: CcdImage.h:93
Polynomial transformation class.
Definition: Gtransfo.h:193
GtransfoLin normalizeCoordinatesTransfo(const Frame &frame)
Returns the transformation that maps the input frame along both axes to [-1,1].
Definition: Gtransfo.cc:724
STL class.
Gtransfo const & getTransfo(CcdImage const &ccdImage) const
Access to mappings.
unsigned assignIndices(unsigned firstIndex, std::string const &whatToFit)
Positions the various parameter sets into the parameter vector, starting at firstIndex.
unsigned getNpar() const
Number of parameters in total.
rectangle with sides parallel to axes.
Definition: Frame.h:19
Class for a simple mapping implementing a generic Gtransfo.
SimplePolyModel(CcdImageList const &ccdImageList, ProjectionHandler const *projectionHandler, bool initFromWCS, unsigned nNotFit=0, unsigned degree=3)
Sky2TP is just a name, it can be anything.
Point getTangentPoint() const
tangent point coordinates (in degrees)
Definition: Gtransfo.cc:1459
This is a virtual class that allows a lot of freedom in the choice of the projection from "Sky" (wher...
void offsetParams(Eigen::VectorXd const &delta)
Offset the parameters by the provided amounts.
const Mapping * getMapping(CcdImage const &) const
Mapping associated to a given CcdImage.
T move(T... args)
T count(T... args)
unsigned getDegree() const
returns degree
Definition: Gtransfo.h:216
Point getCenter() const
Center of the frame.
Definition: Frame.h:41
T find(T... args)
A do-nothing transformation. It anyway has dummy routines to mimick a Gtransfo.
Definition: Gtransfo.h:150
T size(T... args)
#define LSST_EXCEPT(type,...)
virtual GtransfoLin linearApproximation(const Point &where, const double step=0.01) const
linear (local) approximation.
Definition: Gtransfo.cc:92
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:517
T cbegin(T... args)
GtransfoLin getLinPart() const
The Linear part (corresponding to CD&#39;s and CRPIX&#39;s)
Definition: Gtransfo.cc:1461
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:65
a virtual (interface) class for geometric transformations.
Definition: Gtransfo.h:39
int id
void setDegree(const unsigned degree)
Definition: Gtransfo.cc:464
std::shared_ptr< TanSipPix2RaDec > produceSipWcs(CcdImage const &ccdImage) const
Cook up a SIP WCS.
Handler of an actual image from a single CCD.
Definition: CcdImage.h:35
#define LOG_GET(logger)
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:147
#define LOGLS_ERROR(logger, message)
#define LOGLS_WARN(logger, message)
int getNpar() const
total number of parameters
Definition: Gtransfo.h:219