LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
multiShapeletFunction.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 #include "pybind11/pybind11.h"
23 #include "pybind11/stl.h"
24 
25 #include "ndarray/pybind11.h"
26 
28 
29 namespace py = pybind11;
30 using namespace pybind11::literals;
31 
32 namespace lsst {
33 namespace shapelet {
34 
35 namespace {
36 
37 template <typename PyClass>
38 void declareMultiShapeletFunctionMembers(PyClass &cls) {
39  using Class = MultiShapeletFunction;
40 
41  cls.def(py::init<>());
42  cls.def(py::init<Class const &>());
43  cls.def(py::init<typename Class::ComponentList const &>());
44  cls.def(py::init<ShapeletFunction const &>());
45 
46  cls.def("getComponents", [](Class const &self) {
47  auto &components = self.getComponents();
48 
49  py::tuple t(components.size());
50  for (size_t i = 0; i < components.size(); ++i) {
51  t[i] = py::cast(components[i]);
52  }
53 
54  return t;
55  });
56  cls.def("addComponent",
57  [](Class &self, typename Class::Component const &c) { self.getComponents().push_back(c); });
58  cls.def("normalize", &Class::normalize, "value"_a = 1.0);
59  cls.def("shiftInPlace", &Class::shiftInPlace);
60  cls.def("transformInPlace", &Class::transformInPlace);
61  cls.def("convolve", (Class (Class::*)(ShapeletFunction const &) const) & Class::convolve);
62  cls.def("convolve", (Class (Class::*)(Class const &) const) & Class::convolve);
63  cls.def("evaluate", &Class::evaluate);
64 }
65 
66 template <typename PyClass>
67 void declareMultiShapeletFunctionEvaluatorMembers(PyClass &cls) {
68  using Class = MultiShapeletFunctionEvaluator;
69 
70  cls.def(py::init<MultiShapeletFunction const &>());
71 
72  cls.def("__call__", (double (Class::*)(double, double) const) & Class::operator());
73  cls.def("__call__", (double (Class::*)(geom::Point2D const &) const) & Class::operator());
74  cls.def("__call__", (double (Class::*)(geom::Extent2D const &) const) & Class::operator());
75  cls.def("__call__",
76  (ndarray::Array<double, 1, 1> (Class::*)(ndarray::Array<double const, 1> const &,
77  ndarray::Array<double const, 1> const &) const) &
78  Class::operator());
79 
80  cls.def("addToImage",
81  (void (Class::*)(ndarray::Array<double, 2, 1> const &, geom::Point2I const &) const) &
83  "array"_a, "xy0"_a = geom::Point2I());
84  cls.def("addToImage", (void (Class::*)(afw::image::Image<double> &) const) & Class::addToImage,
85  "image"_a);
86 
87  cls.def("integrate", &Class::integrate);
88  cls.def("computeMoments", &Class::computeMoments);
89  cls.def("update", &Class::update);
90 }
91 
92 } // <anonymous>
93 
94 PYBIND11_MODULE(multiShapeletFunction, mod) {
95  py::module::import("lsst.afw.geom");
96  py::module::import("lsst.afw.image");
97 
98  py::class_<MultiShapeletFunction, std::shared_ptr<MultiShapeletFunction>> clsMultiShapeletFunction(
99  mod, "MultiShapeletFunction");
100  py::class_<MultiShapeletFunctionEvaluator, std::shared_ptr<MultiShapeletFunctionEvaluator>>
101  clsMultiShapeletFunctionEvaluator(mod, "MultiShapeletFunctionEvaluator");
102 
103  declareMultiShapeletFunctionMembers(clsMultiShapeletFunction);
104  declareMultiShapeletFunctionEvaluatorMembers(clsMultiShapeletFunctionEvaluator);
105 }
106 
107 } // shapelet
108 } // lsst
table::Key< table::Array< int > > components
py::class_< ProductBoundedField, std::shared_ptr< ProductBoundedField >, BoundedField > PyClass
PYBIND11_MODULE(multiShapeletFunction, mod)
A base class for image defects.
void addToImage(boost::shared_ptr< afw::image::Image< double > > image, std::vector< boost::shared_ptr< afw::image::Image< double > >> const &imgVector, std::vector< double > const &weightVector)
Definition: CoaddPsf.cc:192
UnaryFunctionT::result_type integrate(UnaryFunctionT func, typename UnaryFunctionT::argument_type const a, typename UnaryFunctionT::argument_type const b, double eps=1.0e-6)
The 1D integrator.
Definition: Integrate.h:886