LSST Applications g04e9c324dd+8c5ae1fdc5,g134cb467dc+b203dec576,g18429d2f64+358861cd2c,g199a45376c+0ba108daf9,g1fd858c14a+dd066899e3,g262e1987ae+ebfced1d55,g29ae962dfc+72fd90588e,g2cef7863aa+aef1011c0b,g35bb328faa+8c5ae1fdc5,g3fd5ace14f+b668f15bc5,g4595892280+3897dae354,g47891489e3+abcf9c3559,g4d44eb3520+fb4ddce128,g53246c7159+8c5ae1fdc5,g67b6fd64d1+abcf9c3559,g67fd3c3899+1f72b5a9f7,g74acd417e5+cb6b47f07b,g786e29fd12+668abc6043,g87389fa792+8856018cbb,g89139ef638+abcf9c3559,g8d7436a09f+bcf525d20c,g8ea07a8fe4+9f5ccc88ac,g90f42f885a+6054cc57f1,g97be763408+06f794da49,g9dd6db0277+1f72b5a9f7,ga681d05dcb+7e36ad54cd,gabf8522325+735880ea63,gac2eed3f23+abcf9c3559,gb89ab40317+abcf9c3559,gbf99507273+8c5ae1fdc5,gd8ff7fe66e+1f72b5a9f7,gdab6d2f7ff+cb6b47f07b,gdc713202bf+1f72b5a9f7,gdfd2d52018+8225f2b331,ge365c994fd+375fc21c71,ge410e46f29+abcf9c3559,geaed405ab2+562b3308c0,gf9a733ac38+8c5ae1fdc5,w.2025.35
LSST Data Management Base Package
Loading...
Searching...
No Matches
sdssShape.cc
Go to the documentation of this file.
1/*
2 * LSST Data Management System
3 * Copyright 2008-2017 AURA/LSST.
4 *
5 * This product includes software developed by the
6 * LSST Project (http://www.lsst.org/).
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"
25
26#include <memory>
27
30
34
35namespace py = pybind11;
36using namespace pybind11::literals;
37
38namespace lsst {
39namespace meas {
40namespace base {
41
42namespace {
43
44using PyShapeControl = py::class_<SdssShapeControl>;
45// TODO decide if we need to mention afw::table::FunctorKey<SdssShapeResult>
46// and if so, wrap it; if not, document that here
47using PyShapeResultKey = py::classh<SdssShapeResultKey>;
48using PyShapeResult = py::classh<SdssShapeResult, ShapeResult,
50using PyShapeAlgorithm = py::classh<SdssShapeAlgorithm, SimpleAlgorithm>;
51using PyShapeTransform = py::classh<SdssShapeTransform, BaseTransform>;
52
53PyShapeControl declareShapeControl(lsst::cpputils::python::WrapperCollection &wrappers) {
54 return wrappers.wrapType(PyShapeControl(wrappers.module, "SdssShapeControl"), [](auto &mod, auto &cls) {
55
56 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, background);
57 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, maxIter);
58 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, maxShift);
59 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, tol1);
60 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, tol2);
61 LSST_DECLARE_CONTROL_FIELD(cls, SdssShapeControl, doMeasurePsf);
62
63cls.def(py::init<>());
64 });
65}
66
67void declareShapeResultKey(lsst::cpputils::python::WrapperCollection &wrappers) {
68 wrappers.wrapType(PyShapeResultKey(wrappers.module, "SdssShapeResultKey"), [](auto &mod, auto &cls) {
69 // TODO decide whether to wrap default constructor and do it or document why not
70 cls.def(py::init<afw::table::SubSchema const &>(), "subSchema"_a);
71
72 cls.def_static("addFields", &FluxResultKey::addFields, "schema"_a, "name"_a, "doMeasurePsf"_a);
73
74 cls.def("__eq__", &SdssShapeResultKey::operator==, py::is_operator());
75 cls.def("__ne__", &SdssShapeResultKey::operator!=, py::is_operator());
76
77 cls.def("get", &SdssShapeResultKey::get, "record"_a);
78 cls.def("set", &SdssShapeResultKey::set, "record"_a, "value"_a);
79 cls.def("getPsfShape", &SdssShapeResultKey::getPsfShape, "record"_a);
80 cls.def("setPsfShape", &SdssShapeResultKey::setPsfShape, "record"_a, "value"_a);
81 cls.def("isValid", &SdssShapeResultKey::isValid);
82 cls.def("getFlagHandler", &SdssShapeResultKey::getFlagHandler);
83 });
84}
85
86template <typename ImageT>
87void declareComputeMethods(PyShapeAlgorithm &cls) {
88 cls.def_static("computeAdaptiveMoments",(SdssShapeResult(*)(ImageT const &, geom::Point2D const &, bool, SdssShapeControl const &)) &
90 "image"_a, "position"_a, "negative"_a = false, "ctrl"_a = SdssShapeControl());
91 cls.def_static("computeFixedMomentsFlux",(FluxResult(*)(ImageT const &, afw::geom::ellipses::Quadrupole const &, geom::Point2D const &)) &
93 "image"_a, "shape"_a, "position"_a);
94}
95
96PyShapeAlgorithm declareShapeAlgorithm(lsst::cpputils::python::WrapperCollection &wrappers) {
97 return wrappers.wrapType(PyShapeAlgorithm(wrappers.module, "SdssShapeAlgorithm"), [](auto &mod, auto &cls) {
98 cls.attr("FAILURE") = py::cast(SdssShapeAlgorithm::FAILURE);
99 cls.attr("UNWEIGHTED_BAD") = py::cast(SdssShapeAlgorithm::UNWEIGHTED_BAD);
100 cls.attr("UNWEIGHTED") = py::cast(SdssShapeAlgorithm::UNWEIGHTED);
101 cls.attr("SHIFT") = py::cast(SdssShapeAlgorithm::SHIFT);
102 cls.attr("MAXITER") = py::cast(SdssShapeAlgorithm::MAXITER);
103 cls.attr("PSF_SHAPE_BAD") = py::cast(SdssShapeAlgorithm::PSF_SHAPE_BAD);
104
105 cls.def(py::init<SdssShapeAlgorithm::Control const &, std::string const &, afw::table::Schema &>(),
106 "ctrl"_a, "name"_a, "schema"_a);
107
108 declareComputeMethods<afw::image::Image<int>>(cls);
109 declareComputeMethods<afw::image::Image<float>>(cls);
110 declareComputeMethods<afw::image::Image<double>>(cls);
111 declareComputeMethods<afw::image::MaskedImage<int>>(cls);
112 declareComputeMethods<afw::image::MaskedImage<float>>(cls);
113 declareComputeMethods<afw::image::MaskedImage<double>>(cls);
114
115 cls.def("measure", &SdssShapeAlgorithm::measure, "measRecord"_a, "exposure"_a);
116 cls.def("fail", &SdssShapeAlgorithm::fail, "measRecord"_a, "error"_a = nullptr);
117 });
118}
119
120void declareShapeResult(lsst::cpputils::python::WrapperCollection &wrappers) {
121 wrappers.wrapType(PyShapeResult(wrappers.module, "SdssShapeResult"), [](auto &mod, auto &cls) {
122 cls.def(py::init<>());
123 cls.def_readwrite("instFlux_xx_Cov", &SdssShapeResult::instFlux_xx_Cov);
124 cls.def_readwrite("instFlux_yy_Cov", &SdssShapeResult::instFlux_yy_Cov);
125 cls.def_readwrite("instFlux_xy_Cov", &SdssShapeResult::instFlux_xy_Cov);
126 cls.def_readwrite("flags", &SdssShapeResult::flags);
127
128 // TODO this method says it's a workaround for Swig which doesn't understand std::bitset
129 cls.def("getFlag", (bool (SdssShapeResult::*)(unsigned int) const) &SdssShapeResult::getFlag,
130 "index"_a);
131 cls.def("getFlag", (bool (SdssShapeResult::*)(std::string const &name) const) &SdssShapeResult::getFlag,
132 "name"_a);
133 });
134}
135
136PyShapeTransform declareShapeTransform(lsst::cpputils::python::WrapperCollection &wrappers) {
137 return wrappers.wrapType(PyShapeTransform(wrappers.module, "SdssShapeTransform"), [](auto &mod, auto &cls) {
138 cls.def(py::init<SdssShapeTransform::Control const &, std::string const &, afw::table::SchemaMapper &>(),
139 "ctrl"_a, "name"_a, "mapper"_a);
140
141 cls.def("__call__", &SdssShapeTransform::operator(), "inputCatalog"_a, "outputCatalog"_a, "wcs"_a,
142 "photoCalib"_a);
143 });
144}
145
146} // namespace
147
149 auto clsShapeControl = declareShapeControl(wrappers);
150 declareShapeResultKey(wrappers);
151 auto clsShapeAlgorithm = declareShapeAlgorithm(wrappers);
152 declareShapeResult(wrappers);
153 auto clsShapeTransform = declareShapeTransform(wrappers);
154
155 clsShapeAlgorithm.attr("Control") = clsShapeControl;
156 clsShapeTransform.attr("Control") = clsShapeControl;
157
159 clsShapeAlgorithm, clsShapeControl, clsShapeTransform);
160}
161
162} // namespace base
163} // namespace meas
164} // namespace lsst
A helper class for subdividing pybind11 module across multiple translation units (i....
Definition python.h:242
static Result computeAdaptiveMoments(ImageT const &image, geom::Point2D const &position, bool negative=false, Control const &ctrl=Control())
Compute the adaptive Gaussian-weighted moments of an image.
static FluxResult computeFixedMomentsFlux(ImageT const &image, afw::geom::ellipses::Quadrupole const &shape, geom::Point2D const &position)
Compute the instFlux within a fixed Gaussian aperture.
Definition SdssShape.cc:834
A C++ control class to handle SdssShapeAlgorithm's configuration.
Definition SdssShape.h:52
Result object SdssShapeAlgorithm.
Definition SdssShape.h:224
Point< double, 2 > Point2D
Definition Point.h:324
void declareAlgorithm(PyAlg &clsAlgorithm)
Wrap the implicit API used by meas_base's algorithms.
Definition python.h:87
void wrapSsdsShape(WrapperCollection &)
Definition sdssShape.cc:148
A reusable struct for centroid measurements.
A reusable result struct for instFlux measurements.
A reusable struct for moments-based shape measurements.