lsst.jointcal  16.0-15-g8e16a51+8
SimpleAstrometryModel.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 
4 #include "astshim.h"
5 #include "lsst/log/Log.h"
11 #include "lsst/pex/exceptions.h"
12 #include "lsst/jointcal/Gtransfo.h"
13 
14 namespace {
15 LOG_LOGGER _log = LOG_GET("jointcal.SimpleAstrometryModel");
16 }
17 
18 namespace lsst {
19 namespace jointcal {
20 
23  bool initFromWcs, unsigned nNotFit, unsigned order)
24  : _sky2TP(projectionHandler)
25 
26 {
27  unsigned count = 0;
28 
29  for (auto i = ccdImageList.cbegin(); i != ccdImageList.cend(); ++i, ++count) {
30  const CcdImage &im = **i;
31  if (count < nNotFit) {
33  id->setIndex(-1); // non sense, because it has no parameters
34  _myMap[im.getHashKey()] = std::move(id);
35  } else
36  // Given how AssignIndices works, only the SimplePolyMappings
37  // will actually be fitted, as nNotFit requests.
38  {
39  // first check that there are enough measurements for the requested polynomial order.
40  size_t nObj = im.getCatalogForFit().size();
41  if (nObj == 0) {
42  LOGLS_WARN(_log, "Empty catalog from image: " << im.getName());
43  continue;
44  }
45  GtransfoPoly pol(order);
46  if (pol.getOrder() > 0) // if not, it cannot be decreased
47  {
48  while (unsigned(pol.getNpar()) > 2 * nObj) {
49  LOGLS_WARN(_log, "Reducing polynomial order from "
50  << pol.getOrder() << ", due to too few sources (" << nObj
51  << " vs. " << pol.getNpar() << " parameters)");
52  pol.setOrder(pol.getOrder() - 1);
53  }
54  }
55  /* We have to center and normalize the coordinates so that
56  the fit matrix is not too ill-conditionned. Basically, x
57  and y in pixels are mapped to [-1,1]. When the
58  transformation of SimplePolyMapping transformation is
59  accessed, the combination of the normalization and the
60  fitted transformation is returned, so that the trick
61  remains hidden
62  */
63  const Frame &frame = im.getImageFrame();
64  GtransfoLin shiftAndNormalize = normalizeCoordinatesTransfo(frame);
65  if (initFromWcs) {
66  pol = GtransfoPoly(im.getPix2TangentPlane(), frame, order);
67  pol = pol * shiftAndNormalize.inverted();
68  }
69  _myMap[im.getHashKey()] =
70  std::unique_ptr<SimpleGtransfoMapping>(new SimplePolyMapping(shiftAndNormalize, pol));
71  }
72  }
73 }
74 
76  return findMapping(ccdImage);
77 }
78 
79 unsigned SimpleAstrometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
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 SimpleAstrometryModel::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) i.second->freezeErrorTransform();
103 }
104 
106  int total = 0;
107  for (auto &i : _myMap) {
108  total += i.second->getNpar();
109  }
110  return total;
111 }
112 
114  return dynamic_cast<const SimplePolyMapping *>(findMapping(ccdImage))->getTransfo();
115 }
116 
118  auto proj = std::dynamic_pointer_cast<const TanRaDec2Pix>(getSky2TP(ccdImage));
119  jointcal::Point tangentPoint(proj->getTangentPoint());
120 
121  auto polyMap = getTransfo(ccdImage).toAstMap(ccdImage.getImageFrame());
122  ast::Frame pixelFrame(2, "Domain=PIXELS");
123  ast::Frame iwcFrame(2, "Domain=IWC");
124 
125  // make a basic SkyWcs and extract the IWC portion
126  auto iwcToSkyWcs = afw::geom::makeSkyWcs(
127  afw::geom::Point2D(0, 0),
128  afw::geom::SpherePoint(tangentPoint.x, tangentPoint.y, afw::geom::degrees),
129  afw::geom::makeCdMatrix(1.0 * afw::geom::degrees, 0 * afw::geom::degrees, true));
130  auto iwcToSkyMap = iwcToSkyWcs->getFrameDict()->getMapping("PIXELS", "SKY");
131  auto skyFrame = iwcToSkyWcs->getFrameDict()->getFrame("SKY");
132 
133  ast::FrameDict frameDict(pixelFrame, *polyMap, iwcFrame);
134  frameDict.addFrame("IWC", *iwcToSkyMap, *skyFrame);
135  return std::make_shared<afw::geom::SkyWcs>(frameDict);
136 }
137 
138 AstrometryMapping *SimpleAstrometryModel::findMapping(CcdImage const &ccdImage) const {
139  auto i = _myMap.find(ccdImage.getHashKey());
140  if (i == _myMap.end())
142  "SimpleAstrometryModel cannot find CcdImage " + ccdImage.getName());
143  return i->second.get();
144 }
145 
146 } // namespace jointcal
147 } // namespace lsst
#define LOGLS_WARN(logger, message)
implements the linear transformations (6 real coefficients).
Definition: Gtransfo.h:391
GtransfoLin inverted() const
returns the inverse: T1 = T2.inverted();
Definition: Gtransfo.cc:1178
Gtransfo const & getTransfo(CcdImage const &ccdImage) const
Access to mappings.
SimpleAstrometryModel(CcdImageList const &ccdImageList, const std::shared_ptr< ProjectionHandler const > projectionHandler, bool initFromWCS, unsigned nNotFit=0, unsigned order=3)
Sky2TP is just a name, it can be anything.
A point in a plane.
Definition: Point.h:13
const AstrometryMapping * getMapping(CcdImage const &) const override
Mapping associated to a given CcdImage.
void offsetParams(Eigen::VectorXd const &delta) override
Offset the parameters by the provided amounts (by -delta).
void setOrder(const unsigned order)
Sets the polynomial order (the highest sum of exponents of the largest monomial). ...
Definition: Gtransfo.cc:501
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:56
Mapping implementation for a polynomial transformation.
T cend(T... args)
Gtransfo const * getPix2TangentPlane() const
Definition: CcdImage.h:115
Polynomial transformation class.
Definition: Gtransfo.h:253
GtransfoLin normalizeCoordinatesTransfo(const Frame &frame)
Returns the transformation that maps the input frame along both axes to [-1,1].
Definition: Gtransfo.cc:764
STL class.
std::shared_ptr< afw::geom::SkyWcs > makeSkyWcs(CcdImage const &ccdImage) const override
Make a SkyWcs that contains this model.
unsigned getNpar() const
Number of parameters in total.
const std::shared_ptr< Gtransfo const > getSky2TP(CcdImage const &ccdImage) const override
the mapping of sky coordinates (i.e.
rectangle with sides parallel to axes.
Definition: Frame.h:15
Class for a simple mapping implementing a generic Gtransfo.
void addFrame(int iframe, Mapping const &map, Frame const &frame) override
int getTotalParameters() const override
Return the total number of parameters in this model.
virtual std::shared_ptr< ast::Mapping > toAstMap(jointcal::Frame const &domain) const
Create an equivalent AST mapping for this transformation, including an analytic inverse if possible...
Definition: Gtransfo.h:166
T dynamic_pointer_cast(T... args)
T move(T... args)
T count(T... args)
unsigned getOrder() const
Returns the polynomial order.
Definition: Gtransfo.h:279
T find(T... args)
A do-nothing transformation. It anyway has dummy routines to mimick a Gtransfo.
Definition: Gtransfo.h:196
T size(T... args)
#define LSST_EXCEPT(type,...)
STL class.
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
Definition: Gtransfo.h:653
T cbegin(T... args)
MeasuredStarList const & getCatalogForFit() const
Gets the catalog to be used for fitting, which may have been cleaned-up.
Definition: CcdImage.h:71
a virtual (interface) class for geometric transformations.
Definition: Gtransfo.h:42
std::shared_ptr< SkyWcs > makeSkyWcs(daf::base::PropertySet &metadata, bool strip=false)
int id
Eigen::Matrix2d makeCdMatrix(lsst::geom::Angle const &scale, lsst::geom::Angle const &orientation=0 *lsst::geom::degrees, bool flipX=false)
virtual class needed in the abstraction of the distortion model
CcdImageKey getHashKey() const
Definition: CcdImage.h:128
Handler of an actual image from a single CCD.
Definition: CcdImage.h:41
#define LOG_GET(logger)
unsigned assignIndices(std::string const &whatToFit, unsigned firstIndex) override
Positions the various parameter sets into the parameter vector, starting at firstIndex.
#define LOGLS_ERROR(logger, message)
int getNpar() const override
total number of parameters
Definition: Gtransfo.h:293
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:166