LSST Applications g0f08755f38+f106ab12ff,g12f32b3c4e+379379c50c,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g1a0ca8cf93+eaabd3c39c,g28da252d5a+0daa941822,g2bbee38e9b+c076ca1133,g2bc492864f+c076ca1133,g3156d2b45e+41e33cbcdc,g347aa1857d+c076ca1133,g35bb328faa+a8ce1bb630,g3a166c0a6a+c076ca1133,g3e281a1b8c+b162652f75,g414038480c+6cfc39b08c,g41af890bb2+c487125eda,g5fbc88fb19+17cd334064,g6b1c1869cb+eff3e9fdae,g781aacb6e4+a8ce1bb630,g80478fca09+a0d7f63753,g82479be7b0+052d678814,g858d7b2824+f106ab12ff,g89c8672015+ff349045bf,g9125e01d80+a8ce1bb630,g9726552aa6+3a8748fa0c,ga5288a1d22+d836312400,gb58c049af0+d64f4d3760,gc28159a63d+c076ca1133,gcf0d15dbbd+3012552f84,gd7358e8bfb+8e90b07072,gda3e153d99+f106ab12ff,gda6a2b7d83+3012552f84,gdaeeff99f8+1711a396fd,ge2409df99d+e4ed96189f,ge79ae78c31+c076ca1133,gf0baf85859+50e8a91c63,gf3967379c6+f3639e9197,gfb92a5be7c+f106ab12ff,gfec2e1e490+bdc87d04ab,w.2024.47
LSST Data Management Base Package
Loading...
Searching...
No Matches
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"
36
37namespace py = pybind11;
38using namespace pybind11::literals;
39
40namespace lsst {
41namespace meas {
42namespace extensions {
43namespace photometryKron {
44
45namespace {
46
47void 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
64void 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
91using PyKronAperture = py::class_<KronAperture>;
92
99template <typename ImageT>
100void 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
107void 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
132PYBIND11_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
std::shared_ptr< KronAperture > transform(geom::AffineTransform const &trans) const
Transform a Kron Aperture to a different frame.
static afw::geom::ellipses::Axes getKronAxes(afw::geom::ellipses::Axes const &shape, geom::LinearTransform const &transformation, double const radius)
Determine Kron axes from a reference image.
static meas::base::FlagDefinition const BAD_RADIUS
static meas::base::FlagDefinition const USED_MINIMUM_RADIUS
static meas::base::FlagDefinition const FAILURE
static meas::base::FlagDefinition const SMALL_RADIUS
KronFluxControl Control
A typedef to the Control object for this algorithm, defined above.
virtual void fail(afw::table::SourceRecord &measRecord, meas::base::MeasurementError *error=NULL) const
Handle an exception thrown by the current algorithm by setting flags in the given record.
virtual void measure(afw::table::SourceRecord &measRecord, afw::image::Exposure< float > const &exposure) const
Called to measure a single child source in an image.
static meas::base::FlagDefinition const BAD_SHAPE
static meas::base::FlagDefinitionList const & getFlagDefinitions()
virtual void measureForced(afw::table::SourceRecord &measRecord, afw::image::Exposure< float > const &exposure, afw::table::SourceRecord const &refRecord, afw::geom::SkyWcs const &refWcs) const
Called to measure a single child source in an image.
static meas::base::FlagDefinition const EDGE
static meas::base::FlagDefinition const NO_MINIMUM_RADIUS
static meas::base::FlagDefinition const NO_FALLBACK_RADIUS
static meas::base::FlagDefinition const BAD_SHAPE_NO_PSF
static meas::base::FlagDefinition const USED_PSF_RADIUS
T nanf(T... args)
#define LSST_DECLARE_CONTROL_FIELD(WRAPPER, CLASS, NAME)
Macro used to wrap fields declared by LSST_CONTROL_FIELD using Pybind11.
Definition python.h:50