LSSTApplications  17.0+124,17.0+14,17.0+73,18.0.0+37,18.0.0+80,18.0.0-4-g68ffd23+4,18.1.0-1-g0001055+12,18.1.0-1-g03d53ef+5,18.1.0-1-g1349e88+55,18.1.0-1-g2505f39+44,18.1.0-1-g5315e5e+4,18.1.0-1-g5e4b7ea+14,18.1.0-1-g7e8fceb+4,18.1.0-1-g85f8cd4+48,18.1.0-1-g8ff0b9f+4,18.1.0-1-ga2c679d+1,18.1.0-1-gd55f500+35,18.1.0-10-gb58edde+2,18.1.0-11-g0997b02+4,18.1.0-13-gfe4edf0b+12,18.1.0-14-g259bd21+21,18.1.0-19-gdb69f3f+2,18.1.0-2-g5f9922c+24,18.1.0-2-gd3b74e5+11,18.1.0-2-gfbf3545+32,18.1.0-26-g728bddb4+5,18.1.0-27-g6ff7ca9+2,18.1.0-3-g52aa583+25,18.1.0-3-g8ea57af+9,18.1.0-3-gb69f684+42,18.1.0-3-gfcaddf3+6,18.1.0-32-gd8786685a,18.1.0-4-gf3f9b77+6,18.1.0-5-g1dd662b+2,18.1.0-5-g6dbcb01+41,18.1.0-6-gae77429+3,18.1.0-7-g9d75d83+9,18.1.0-7-gae09a6d+30,18.1.0-9-gc381ef5+4,w.2019.45
LSSTDataManagementBasePackage
photometryKron.cc
Go to the documentation of this file.
1 /*
2  * LSST Data Management System
3  *
4  * This product includes software developed by the
5  * LSST Project (http://www.lsst.org/).
6  * See the COPYRIGHT file
7  *
8  * This program is free software: you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation, either version 3 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the LSST License Statement and
19  * the GNU General Public License along with this program. If not,
20  * see <https://www.lsstcorp.org/LegalNotices/>.
21  */
22 
23 #include "pybind11/pybind11.h"
24 #include "pybind11/stl.h"
25 
26 #include <cstdint>
27 
28 #include "lsst/pex/config/python.h" // defines LSST_DECLARE_CONTROL_FIELD
29 #include "lsst/geom/Point.h"
30 #include "lsst/afw/geom/ellipses.h"
33 #include "lsst/afw/table/Source.h"
36 
37 namespace py = pybind11;
38 using namespace pybind11::literals;
39 
40 namespace lsst {
41 namespace meas {
42 namespace extensions {
43 namespace photometryKron {
44 
45 namespace {
46 
47 void declareKronFluxControl(py::module &mod) {
48  py::class_<KronFluxControl> cls(mod, "KronFluxControl");
49 
50  cls.def(py::init<>());
51 
52  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, fixed);
53  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, nSigmaForRadius);
54  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, nIterForRadius);
55  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, nRadiusForFlux);
56  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, maxSincRadius);
57  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, minimumRadius);
58  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, enforceMinimumRadius);
59  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, useFootprintRadius);
60  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, smoothingSigma);
61  LSST_DECLARE_CONTROL_FIELD(cls, KronFluxControl, refRadiusName);
62 }
63 
64 void declareKronFluxAlgorithm(py::module &mod) {
65  py::class_<KronFluxAlgorithm, std::shared_ptr<KronFluxAlgorithm>, base::SimpleAlgorithm> cls(
66  mod, "KronFluxAlgorithm");
67 
68  cls.def_static("getFlagDefinitions", &KronFluxAlgorithm::getFlagDefinitions,
69  py::return_value_policy::copy);
70  cls.attr("FAILURE") = py::cast(KronFluxAlgorithm::FAILURE);
71  cls.attr("EDGE") = py::cast(KronFluxAlgorithm::EDGE);
72  cls.attr("BAD_SHAPE_NO_PSF") = py::cast(KronFluxAlgorithm::BAD_SHAPE_NO_PSF);
73  cls.attr("NO_MINIMUM_RADIUS") = py::cast(KronFluxAlgorithm::NO_MINIMUM_RADIUS);
74  cls.attr("NO_FALLBACK_RADIUS") = py::cast(KronFluxAlgorithm::NO_FALLBACK_RADIUS);
75  cls.attr("BAD_RADIUS") = py::cast(KronFluxAlgorithm::BAD_RADIUS);
76  cls.attr("USED_MINIMUM_RADIUS") = py::cast(KronFluxAlgorithm::USED_MINIMUM_RADIUS);
77  cls.attr("USED_PSF_RADIUS") = py::cast(KronFluxAlgorithm::USED_PSF_RADIUS);
78  cls.attr("SMALL_RADIUS") = py::cast(KronFluxAlgorithm::SMALL_RADIUS);
79  cls.attr("BAD_SHAPE") = py::cast(KronFluxAlgorithm::BAD_SHAPE);
80 
81  cls.def(py::init<KronFluxAlgorithm::Control const &, std::string const &, afw::table::Schema &,
82  daf::base::PropertySet &>(),
83  "ctrl"_a, "name"_a, "schema"_a, "metadata"_a);
84 
85  cls.def("measure", &KronFluxAlgorithm::measure, "measRecord"_a, "exposure"_a);
86  cls.def("measureForced", &KronFluxAlgorithm::measureForced, "measRecord"_a, "exposure"_a, "refRecord"_a,
87  "refWcs"_a);
88  cls.def("fail", &KronFluxAlgorithm::fail, "measRecord"_a, "error"_a = NULL);
89 }
90 
91 using PyKronAperture = py::class_<KronAperture>;
92 
99 template <typename ImageT>
100 void declareKronApertureTemplatedMethods(PyKronAperture &cls) {
101  cls.def_static("determineRadius", &KronAperture::determineRadius<ImageT>, "image"_a, "axes"_a, "center"_a,
102  "ctrl"_a);
103  cls.def("measureFlux", &KronAperture::measureFlux<ImageT>, "image"_a, "nRadiusForFlux"_a,
104  "maxSincRadius"_a);
105 }
106 
107 void declareKronAperture(py::module &mod) {
108  PyKronAperture cls(mod, "KronAperture");
109 
110  cls.def(py::init<geom::Point2D const &, afw::geom::ellipses::BaseCore const &, float>(), "center"_a,
111  "core"_a, "radiusForRadius"_a = std::nanf(""));
112  cls.def(py::init<afw::table::SourceRecord const &, float>(), "source"_a,
113  "radiusForRadius"_a = std::nanf(""));
114  cls.def(py::init<afw::table::SourceRecord const &, geom::AffineTransform const &, double, float>(),
115  "reference"_a, "refToMeas"_a, "radius"_a, "radiusForRadius"_a = std::nanf(""));
116 
117  cls.def_static("getKronAxes", &KronAperture::getKronAxes, "shape"_a, "transformation"_a, "radius"_a);
118 
119  cls.def("getX", &KronAperture::getX);
120  cls.def("getY", &KronAperture::getY);
121  cls.def("getRadiusForRadius", &KronAperture::getRadiusForRadius);
122  cls.def("getCenter", &KronAperture::getCenter);
123  cls.def("getAxes", (afw::geom::ellipses::Axes & (KronAperture::*)()) & KronAperture::getAxes,
124  py::return_value_policy::reference_internal);
125  cls.def("transform", &KronAperture::transform, "trans"_a);
126 
127  declareKronApertureTemplatedMethods<afw::image::MaskedImage<float>>(cls);
128 }
129 
130 } // <anonymous>
131 
132 PYBIND11_MODULE(photometryKron, mod) {
133  py::module::import("lsst.afw.geom");
134  py::module::import("lsst.afw.image");
135  py::module::import("lsst.afw.table");
136  py::module::import("lsst.daf.base");
137 
138  declareKronFluxControl(mod);
139  declareKronFluxAlgorithm(mod);
140  declareKronAperture(mod);
141 }
142 
143 } // photometryKron
144 } // extensions
145 } // meas
146 } // lsst
def init()
Definition: tests.py:64
STL class.
A base class for image defects.
#define LSST_DECLARE_CONTROL_FIELD(WRAPPER, CLASS, NAME)
Macro used to wrap fields declared by LSST_CONTROL_FIELD using Pybind11.
Definition: python.h:50
T nanf(T... args)
def measure(mi, x, y, size, statistic, stats)
Definition: fringe.py:473
table::Key< int > transform