lsst.jointcal  17.0.1-12-g5e92ed9+4
SimpleAstrometryModel.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 <memory>
26 #include <string>
27 
28 #include "astshim.h"
29 #include "lsst/log/Log.h"
33 #include "lsst/jointcal/CcdImage.h"
35 #include "lsst/pex/exceptions.h"
37 
38 namespace lsst {
39 namespace jointcal {
40 
43  bool initFromWcs, unsigned nNotFit, unsigned order)
44  : AstrometryModel(LOG_GET("jointcal.SimpleAstrometryModel")),
45  _skyToTangentPlane(projectionHandler)
46 
47 {
48  unsigned count = 0;
49 
50  for (auto i = ccdImageList.cbegin(); i != ccdImageList.cend(); ++i, ++count) {
51  const CcdImage &im = **i;
52  if (count < nNotFit) {
55  id->setIndex(-1); // non sense, because it has no parameters
56  _myMap[im.getHashKey()] = std::move(id);
57  } else
58  // Given how AssignIndices works, only the SimplePolyMappings
59  // will actually be fitted, as nNotFit requests.
60  {
61  // first check that there are enough measurements for the requested polynomial order.
62  size_t nObj = im.getCatalogForFit().size();
63  if (nObj == 0) {
64  LOGLS_WARN(_log, "Empty catalog from image: " << im.getName());
65  continue;
66  }
68  if (pol.getOrder() > 0) // if not, it cannot be decreased
69  {
70  while (unsigned(pol.getNpar()) > 2 * nObj) {
71  LOGLS_WARN(_log, "Reducing polynomial order from "
72  << pol.getOrder() << ", due to too few sources (" << nObj
73  << " vs. " << pol.getNpar() << " parameters)");
74  pol.setOrder(pol.getOrder() - 1);
75  }
76  }
77  /* We have to center and normalize the coordinates so that
78  the fit matrix is not too ill-conditionned. Basically, x
79  and y in pixels are mapped to [-1,1]. When the
80  transformation of SimplePolyMapping transformation is
81  accessed, the combination of the normalization and the
82  fitted transformation is returned, so that the trick
83  remains hidden
84  */
85  const Frame &frame = im.getImageFrame();
87  if (initFromWcs) {
88  pol = AstrometryTransformPolynomial(im.getPixelToTangentPlane().get(), frame, order);
89  pol = pol * shiftAndNormalize.inverted();
90  }
91  _myMap[im.getHashKey()] =
93  }
94  }
95 }
96 
98  return findMapping(ccdImage);
99 }
100 
101 unsigned SimpleAstrometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
102  if (whatToFit.find("Distortions") == std::string::npos) {
103  LOGLS_ERROR(_log, "AssignIndices was called and Distortions is *not* in whatToFit.");
104  return 0;
105  }
106  unsigned index = firstIndex;
107  for (auto i = _myMap.begin(); i != _myMap.end(); ++i) {
108  SimplePolyMapping *p = dynamic_cast<SimplePolyMapping *>(&*(i->second));
109  if (!p) continue; // it should be AstrometryTransformIdentity
110  p->setIndex(index);
111  index += p->getNpar();
112  }
113  return index;
114 }
115 
116 void SimpleAstrometryModel::offsetParams(Eigen::VectorXd const &delta) {
117  for (auto &i : _myMap) {
118  auto mapping = i.second.get();
119  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
120  }
121 }
122 
124  for (auto &i : _myMap) i.second->freezeErrorTransform();
125 }
126 
128  int total = 0;
129  for (auto &i : _myMap) {
130  total += i.second->getNpar();
131  }
132  return total;
133 }
134 
136  return dynamic_cast<const SimplePolyMapping *>(findMapping(ccdImage))->getTransform();
137 }
138 
140  auto proj = std::dynamic_pointer_cast<const TanRaDecToPixel>(getSkyToTangentPlane(ccdImage));
141  jointcal::Point tangentPoint(proj->getTangentPoint());
142 
143  auto polyMap = getTransform(ccdImage).toAstMap(ccdImage.getImageFrame());
144  ast::Frame pixelFrame(2, "Domain=PIXELS");
145  ast::Frame iwcFrame(2, "Domain=IWC");
146 
147  // make a basic SkyWcs and extract the IWC portion
148  auto iwcToSkyWcs = afw::geom::makeSkyWcs(
149  afw::geom::Point2D(0, 0),
150  afw::geom::SpherePoint(tangentPoint.x, tangentPoint.y, afw::geom::degrees),
151  afw::geom::makeCdMatrix(1.0 * afw::geom::degrees, 0 * afw::geom::degrees, true));
152  auto iwcToSkyMap = iwcToSkyWcs->getFrameDict()->getMapping("PIXELS", "SKY");
153  auto skyFrame = iwcToSkyWcs->getFrameDict()->getFrame("SKY");
154 
155  ast::FrameDict frameDict(pixelFrame, *polyMap, iwcFrame);
156  frameDict.addFrame("IWC", *iwcToSkyMap, *skyFrame);
157  return std::make_shared<afw::geom::SkyWcs>(frameDict);
158 }
159 
160 AstrometryMapping *SimpleAstrometryModel::findMapping(CcdImage const &ccdImage) const {
161  auto i = _myMap.find(ccdImage.getHashKey());
162  if (i == _myMap.end())
164  "SimpleAstrometryModel cannot find CcdImage " + ccdImage.getName());
165  return i->second.get();
166 }
167 
168 } // namespace jointcal
169 } // namespace lsst
AstrometryTransform const & getTransform(CcdImage const &ccdImage) const
Access to mappings.
void setOrder(const unsigned order)
Sets the polynomial order (the highest sum of exponents of the largest monomial). ...
implements the linear transformations (6 real coefficients).
SimpleAstrometryModel(CcdImageList const &ccdImageList, const std::shared_ptr< ProjectionHandler const > projectionHandler, bool initFromWCS, unsigned nNotFit=0, unsigned order=3)
A point in a plane.
Definition: Point.h:36
LOG_LOGGER _log
lsst.logging instance, to be created by a subclass so that messages have consistent name...
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).
std::shared_ptr< AstrometryTransform > const getPixelToTangentPlane() const
Definition: CcdImage.h:140
std::string getName() const
Return the _name that identifies this ccdImage.
Definition: CcdImage.h:79
Mapping implementation for a polynomial transformation.
unsigned getNpar() const
Number of parameters in total.
Interface between AstrometryFit and the combinations of Mappings from pixels to some tangent plane (a...
const std::shared_ptr< AstrometryTransform const > getSkyToTangentPlane(CcdImage const &ccdImage) const override
the mapping of sky coordinates (i.e.
T cend(T... args)
STL class.
std::shared_ptr< afw::geom::SkyWcs > makeSkyWcs(CcdImage const &ccdImage) const override
Make a SkyWcs that contains this model.
void addFrame(int iframe, Mapping const &map, Frame const &frame) override
rectangle with sides parallel to axes.
Definition: Frame.h:38
AstrometryTransformLinear normalizeCoordinatesTransform(const Frame &frame)
Returns the transformation that maps the input frame along both axes to [-1,1].
Class for a simple mapping implementing a generic AstrometryTransform.
int getTotalParameters() const override
Return the total number of parameters in this model.
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane) ...
T dynamic_pointer_cast(T... args)
table::Key< int > id
T move(T... args)
T count(T... args)
int getNpar() const override
total number of parameters
A do-nothing transformation. It anyway has dummy routines to mimick a AstrometryTransform.
T find(T... args)
T size(T... args)
#define LSST_EXCEPT(type,...)
STL class.
AstrometryTransformLinear inverted() const
returns the inverse: T1 = T2.inverted();
a virtual (interface) class for geometric transformations.
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:94
std::shared_ptr< SkyWcs > makeSkyWcs(daf::base::PropertySet &metadata, bool strip=false)
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:152
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
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...
unsigned assignIndices(std::string const &whatToFit, unsigned firstIndex) override
Positions the various parameter sets into the parameter vector, starting at firstIndex.
unsigned getOrder() const
Returns the polynomial order.
#define LOG_GET(logger)
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:191
#define LOGLS_ERROR(logger, message)
#define LOGLS_WARN(logger, message)