LSSTApplications  16.0-11-g09ed895+2,16.0-11-g12e47bd,16.0-11-g9bb73b2+6,16.0-12-g5c924a4+6,16.0-14-g9a974b3+1,16.0-15-g1417920+1,16.0-15-gdd5ca33+1,16.0-16-gf0259e2,16.0-17-g31abd91+7,16.0-17-g7d7456e+7,16.0-17-ga3d2e9f+13,16.0-18-ga4d4bcb+1,16.0-18-gd06566c+1,16.0-2-g0febb12+21,16.0-2-g9d5294e+69,16.0-2-ga8830df+6,16.0-20-g21842373+7,16.0-24-g3eae5ec,16.0-28-gfc9ea6c+4,16.0-29-ge8801f9,16.0-3-ge00e371+34,16.0-4-g18f3627+13,16.0-4-g5f3a788+20,16.0-4-ga3eb747+10,16.0-4-gabf74b7+29,16.0-4-gb13d127+6,16.0-49-g42e581f7+6,16.0-5-g27fb78a+7,16.0-5-g6a53317+34,16.0-5-gb3f8a4b+87,16.0-6-g9321be7+4,16.0-6-gcbc7b31+42,16.0-6-gf49912c+29,16.0-7-gd2eeba5+51,16.0-71-ge89f8615e,16.0-8-g21fd5fe+29,16.0-8-g3a9f023+20,16.0-8-g4734f7a+1,16.0-8-g5858431+3,16.0-9-gf5c1f43+8,master-gd73dc1d098+1,w.2019.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::*)(afw::geom::Point2D const &) const) & Class::operator());
74  cls.def("__call__", (double (Class::*)(afw::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 &, afw::geom::Point2I const &) const) &
83  "array"_a, "xy0"_a = afw::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
PYBIND11_MODULE(multiShapeletFunction, mod)
Point< double, 2 > Point2D
Definition: Point.h:324
A base class for image defects.
Definition: cameraGeom.dox:3
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, ConvolutionControl const &convolutionControl=ConvolutionControl())
Convolve an Image or MaskedImage with a Kernel, setting pixels of an existing output image...
Point< int, 2 > Point2I
Definition: Point.h:321
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
Extent< double, 2 > Extent2D
Definition: Extent.h:400
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