lsst.jointcal  16.0-28-gfc9ea6c+15
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 {
39 LOG_LOGGER _log = LOG_GET("jointcal.SimpleAstrometryModel");
40 }
41 
42 namespace lsst {
43 namespace jointcal {
44 
47  bool initFromWcs, unsigned nNotFit, unsigned order)
48  : _skyToTangentPlane(projectionHandler)
49 
50 {
51  unsigned count = 0;
52 
53  for (auto i = ccdImageList.cbegin(); i != ccdImageList.cend(); ++i, ++count) {
54  const CcdImage &im = **i;
55  if (count < nNotFit) {
58  id->setIndex(-1); // non sense, because it has no parameters
59  _myMap[im.getHashKey()] = std::move(id);
60  } else
61  // Given how AssignIndices works, only the SimplePolyMappings
62  // will actually be fitted, as nNotFit requests.
63  {
64  // first check that there are enough measurements for the requested polynomial order.
65  size_t nObj = im.getCatalogForFit().size();
66  if (nObj == 0) {
67  LOGLS_WARN(_log, "Empty catalog from image: " << im.getName());
68  continue;
69  }
71  if (pol.getOrder() > 0) // if not, it cannot be decreased
72  {
73  while (unsigned(pol.getNpar()) > 2 * nObj) {
74  LOGLS_WARN(_log, "Reducing polynomial order from "
75  << pol.getOrder() << ", due to too few sources (" << nObj
76  << " vs. " << pol.getNpar() << " parameters)");
77  pol.setOrder(pol.getOrder() - 1);
78  }
79  }
80  /* We have to center and normalize the coordinates so that
81  the fit matrix is not too ill-conditionned. Basically, x
82  and y in pixels are mapped to [-1,1]. When the
83  transformation of SimplePolyMapping transformation is
84  accessed, the combination of the normalization and the
85  fitted transformation is returned, so that the trick
86  remains hidden
87  */
88  const Frame &frame = im.getImageFrame();
90  if (initFromWcs) {
91  pol = AstrometryTransformPolynomial(im.getPixelToTangentPlane().get(), frame, order);
92  pol = pol * shiftAndNormalize.inverted();
93  }
94  _myMap[im.getHashKey()] =
96  }
97  }
98 }
99 
101  return findMapping(ccdImage);
102 }
103 
104 unsigned SimpleAstrometryModel::assignIndices(std::string const &whatToFit, unsigned firstIndex) {
105  if (whatToFit.find("Distortions") == std::string::npos) {
106  LOGLS_ERROR(_log, "AssignIndices was called and Distortions is *not* in whatToFit.");
107  return 0;
108  }
109  unsigned index = firstIndex;
110  for (auto i = _myMap.begin(); i != _myMap.end(); ++i) {
111  SimplePolyMapping *p = dynamic_cast<SimplePolyMapping *>(&*(i->second));
112  if (!p) continue; // it should be AstrometryTransformIdentity
113  p->setIndex(index);
114  index += p->getNpar();
115  }
116  return index;
117 }
118 
119 void SimpleAstrometryModel::offsetParams(Eigen::VectorXd const &delta) {
120  for (auto &i : _myMap) {
121  auto mapping = i.second.get();
122  mapping->offsetParams(delta.segment(mapping->getIndex(), mapping->getNpar()));
123  }
124 }
125 
127  for (auto &i : _myMap) i.second->freezeErrorTransform();
128 }
129 
131  int total = 0;
132  for (auto &i : _myMap) {
133  total += i.second->getNpar();
134  }
135  return total;
136 }
137 
139  return dynamic_cast<const SimplePolyMapping *>(findMapping(ccdImage))->getTransform();
140 }
141 
143  auto proj = std::dynamic_pointer_cast<const TanRaDecToPixel>(getSkyToTangentPlane(ccdImage));
144  jointcal::Point tangentPoint(proj->getTangentPoint());
145 
146  auto polyMap = getTransform(ccdImage).toAstMap(ccdImage.getImageFrame());
147  ast::Frame pixelFrame(2, "Domain=PIXELS");
148  ast::Frame iwcFrame(2, "Domain=IWC");
149 
150  // make a basic SkyWcs and extract the IWC portion
151  auto iwcToSkyWcs = afw::geom::makeSkyWcs(
152  afw::geom::Point2D(0, 0),
153  afw::geom::SpherePoint(tangentPoint.x, tangentPoint.y, afw::geom::degrees),
154  afw::geom::makeCdMatrix(1.0 * afw::geom::degrees, 0 * afw::geom::degrees, true));
155  auto iwcToSkyMap = iwcToSkyWcs->getFrameDict()->getMapping("PIXELS", "SKY");
156  auto skyFrame = iwcToSkyWcs->getFrameDict()->getFrame("SKY");
157 
158  ast::FrameDict frameDict(pixelFrame, *polyMap, iwcFrame);
159  frameDict.addFrame("IWC", *iwcToSkyMap, *skyFrame);
160  return std::make_shared<afw::geom::SkyWcs>(frameDict);
161 }
162 
163 AstrometryMapping *SimpleAstrometryModel::findMapping(CcdImage const &ccdImage) const {
164  auto i = _myMap.find(ccdImage.getHashKey());
165  if (i == _myMap.end())
167  "SimpleAstrometryModel cannot find CcdImage " + ccdImage.getName());
168  return i->second.get();
169 }
170 
171 } // namespace jointcal
172 } // namespace lsst
#define LOGLS_WARN(logger, message)
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)
Sky2TP is just a name, it can be anything.
A point in a plane.
Definition: Point.h:36
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.
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.
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.
void addFrame(int iframe, Mapping const &map, Frame const &frame) override
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
#define LOG_GET(logger)
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.
#define LOGLS_ERROR(logger, message)
unsigned getOrder() const
Returns the polynomial order.
Frame const & getImageFrame() const
Frame in pixels.
Definition: CcdImage.h:191