LSSTApplications  17.0+103,17.0+11,17.0+61,18.0.0+13,18.0.0+25,18.0.0+5,18.0.0+52,18.0.0-4-g68ffd23,18.1.0-1-g0001055+8,18.1.0-1-g03d53ef+1,18.1.0-1-g1349e88+28,18.1.0-1-g2505f39+22,18.1.0-1-g380d4d4+27,18.1.0-1-g5315e5e+1,18.1.0-1-g5e4b7ea+10,18.1.0-1-g7e8fceb+1,18.1.0-1-g85f8cd4+23,18.1.0-1-g9a6769a+13,18.1.0-1-ga1a4c1a+22,18.1.0-1-gd55f500+17,18.1.0-12-g42eabe8e+10,18.1.0-14-gd04256d+15,18.1.0-16-g430f6a53+1,18.1.0-17-gd2166b6e4,18.1.0-18-gb5d19ff+1,18.1.0-2-gfbf3545+7,18.1.0-2-gfefb8b5+16,18.1.0-3-g52aa583+13,18.1.0-3-g62b5e86+14,18.1.0-3-g8f4a2b1+17,18.1.0-3-g9bc06b8+7,18.1.0-3-gb69f684+9,18.1.0-4-g1ee41a7+1,18.1.0-5-g6dbcb01+13,18.1.0-5-gc286bb7+3,18.1.0-6-g48bdcd3+2,18.1.0-6-gd05e160+9,18.1.0-7-gc4d902b+2,18.1.0-7-gebc0338+8,18.1.0-9-gae7190a+10,w.2019.38
LSSTDataManagementBasePackage
astrometry_net.cc
Go to the documentation of this file.
1 // -*- lsst-C++ -*-
2 
3 // Astrometry.net include files...
4 extern "C" {
5 #include "astrometry/solver.h"
6 #include "astrometry/index.h"
7 #include "astrometry/multiindex.h"
8 #include "astrometry/starkd.h"
9 #include "astrometry/fitsioutils.h"
10 #include "astrometry/fitstable.h"
11 #include "astrometry/log.h"
12 #include "astrometry/tic.h"
13 #include "astrometry/healpix.h"
14 
15 #undef ATTRIB_FORMAT
16 #undef FALSE
17 #undef TRUE
18 #undef DEG_PER_RAD
19 #undef RAD_PER_DEG
20 #undef logdebug
21 #undef debug
22 }
23 
24 #include <string>
25 #include <sstream>
26 
27 #include "pybind11/pybind11.h"
28 #include "pybind11/stl.h"
29 
30 #include "lsst/log/Log.h"
31 #include "lsst/utils/python.h"
33 
34 namespace py = pybind11;
35 using namespace pybind11::literals;
36 
37 namespace lsst {
38 namespace meas {
39 namespace extensions {
40 namespace astrometryNet {
41 namespace {
42 
43 /*
44  * Wrap index_t
45  *
46  * This class has no constructor; the wrapper is only used for objects returned by MultiIndex::operator[]
47  */
48 static void declareIndex(py::module& mod) {
49  py::class_<index_t> cls(mod, "index_t");
50 
51  cls.def("overlapsScaleRange",
52  [](index_t& self, double qlo, double qhi) { return index_overlaps_scale_range(&self, qlo, qhi); },
53  "qlow"_a, "qhigh"_a);
54  cls.def("reload", [](index_t& self) {
55  if (index_reload(&self)) {
57  os << "Failed to reload multi-index file " << self.indexname;
59  }
60  });
61 
62  cls.def_readonly("indexid", &index_t::indexid);
63  cls.def_readonly("healpix", &index_t::healpix);
64  cls.def_readonly("hpnside", &index_t::hpnside);
65  cls.def_readonly("nstars", &index_t::nstars);
66  cls.def_readonly("nquads", &index_t::nquads);
67 }
68 
72 static void declareMultiIndex(py::module& mod) {
73  py::class_<MultiIndex> cls(mod, "MultiIndex");
74 
75  cls.def(py::init<std::string const&>(), "filepath"_a);
76 
77  cls.def("__getitem__",
78  [](MultiIndex const& self, int i) {
79  auto cind = lsst::utils::python::cppIndex(self.getLength(), i);
80  return self[cind];
81  },
82  py::return_value_policy::reference_internal, py::is_operator());
83 
84  cls.def("addIndex", &MultiIndex::addIndex, "filepath"_a, "metadataOnly"_a);
85  cls.def("isWithinRange", &MultiIndex::isWithinRange, "ra"_a, "dec"_a, "radius"_a);
86  cls.def("unload", &MultiIndex::unload);
87  cls.def_property_readonly("name", &MultiIndex::getName);
88  cls.def("__len__", &MultiIndex::getLength);
89  cls.def("reload", &MultiIndex::reload);
90 }
91 
95 static void declareSolver(py::module& mod) {
96  py::class_<Solver> cls(mod, "Solver");
97 
98  cls.def(py::init<>());
99 
127  cls.def("getCatalog", &Solver::getCatalog, "inds"_a, "ctrCoord"_a, "radius"_a, "idCol"_a,
128  "filterNameList"_a, "magColList"_a, "magErrColList"_a, "starGalCol"_a, "varCol"_a,
129  "uniqueIds"_a = true);
130  cls.def("getSolveStats", &Solver::getSolveStats);
131  cls.def("getWcs", &Solver::getWcs);
132  cls.def("didSolve", &Solver::didSolve);
133  cls.def("run", &Solver::run, "cpulimit"_a);
134  cls.def("getQuadSizeRangeArcsec", &Solver::getQuadSizeRangeArcsec);
135  cls.def("addIndices", &Solver::addIndices, "indices"_a);
136  cls.def("setParity", &Solver::setParity, "setParityFlipped", "parity"_a);
137  cls.def("setMatchThreshold", &Solver::setMatchThreshold, "threshold"_a);
138  cls.def("setPixelScaleRange", &Solver::setPixelScaleRange, "low"_a, "high"_a);
139  cls.def("setRaDecRadius", &Solver::setRaDecRadius, "ra"_a, "dec"_a, "rad"_a);
140  cls.def("setImageSize", &Solver::setImageSize, "width"_a, "height"_a);
141  cls.def("setMaxStars", &Solver::setMaxStars, "maxStars"_a);
142  cls.def("setStars", &Solver::setStars, "sourceCat"_a, "x0"_a, "y0"_a);
143 }
144 
145 // declare logging functions for use by the Python
146 
147 LOG_LOGGER an_log = LOG_GET("meas.astrom.astrometry_net");
148 
149 void an_log_callback(void* baton, enum log_level level, const char* file, int line, const char* func,
150  const char* format, va_list va) {
151  // translate between logging levels
152  int levelmap[5];
153  levelmap[LOG_NONE] = LOG_LVL_FATAL;
154  levelmap[LOG_ERROR] = LOG_LVL_FATAL;
155  levelmap[LOG_MSG] = LOG_LVL_INFO;
156  levelmap[LOG_VERB] = LOG_LVL_DEBUG;
157  levelmap[LOG_ALL] = LOG_LVL_DEBUG;
158  int lsstlevel = levelmap[level];
159 
160  va_list vb;
161  // find out how long the formatted string will be
162  va_copy(vb, va);
163  const int len = vsnprintf(NULL, 0, format, va) + 1; // "+ 1" for the '\0'
164  va_end(va);
165  // allocate a string of the appropriate length
166  char msg[len];
167  (void)vsnprintf(msg, len, format, vb);
168  va_end(vb);
169 
170  // trim trailing \n
171  if (msg[len - 2] == '\n') {
172  msg[len - 2] = '\0';
173  }
174 
175  // not using the LOGS macro because the original location info is wanted
176  if (an_log.isEnabledFor(lsstlevel)) {
177  an_log.logMsg(log4cxx::Level::toLevel(lsstlevel), log4cxx::spi::LocationInfo(file, func, line), msg);
178  }
179 }
180 
182 void start_an_logging() {
183  // NOTE, this has to happen before the log_use_function!
184  log_init(LOG_VERB);
185  log_use_function(an_log_callback, NULL);
186  log_to(NULL);
187 }
188 
190 void finalize() {
191  log_use_function(NULL, NULL);
192  log_to(stdout);
193 }
194 
195 } // namespace
196 
197 PYBIND11_MODULE(astrometry_net, mod) {
198  // code that is run at import time
199  fits_use_error_system();
200  start_an_logging();
201 
202  mod.def("healpixDistance", &healpixDistance, "hp"_a, "nside"_a, "coord"_a);
203 
204  mod.def("an_log_init", [](int level) { log_init(static_cast<log_level>(level)); }, "level"_a);
205 
206  mod.def("an_log_get_level", []() { return static_cast<int>(log_get_level()); });
207  mod.def("an_log_set_level", [](int level) { log_set_level(static_cast<log_level>(level)); }, "level"_a);
208  mod.def("finalize", &finalize);
209 
210  declareMultiIndex(mod);
211  declareIndex(mod);
212  declareSolver(mod);
213 }
214 
215 } // namespace astrometryNet
216 } // namespace extensions
217 } // namespace meas
218 } // namespace lsst
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:167
#define LOG_LOGGER
Definition: Log.h:703
std::size_t cppIndex(std::ptrdiff_t size, std::ptrdiff_t i)
Compute a C++ index from a Python index (negative values count from the end) and range-check.
Definition: python.h:124
LSST DM logging module built on log4cxx.
lsst::geom::Angle healpixDistance(int hp, int nside, lsst::geom::SpherePoint const &coord)
Calculate the distance from coordinates to a healpix.
#define LOG_LVL_INFO
Definition: Log.h:698
A base class for image defects.
#define LOG_LVL_DEBUG
Definition: Log.h:697
T str(T... args)
#define LOG_ERROR(message...)
Log a error-level message to the default logger using a varargs/printf style interface.
Definition: Log.h:315
T vsnprintf(T... args)
def line(points, frame=None, origin=afwImage.PARENT, symbs=False, ctype=None, size=0.5)
Definition: ds9.py:105
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
#define LOG_GET(logger)
Returns a Log object associated with logger.
Definition: Log.h:75
std::ostream * os
Definition: Schema.cc:746
Reports errors that are due to events beyond the control of the program.
Definition: Runtime.h:104
#define LOG_LVL_FATAL
Definition: Log.h:701