LSST Applications  21.0.0+04719a4bac,21.0.0-1-ga51b5d4+f5e6047307,21.0.0-11-g2b59f77+a9c1acf22d,21.0.0-11-ga42c5b2+86977b0b17,21.0.0-12-gf4ce030+76814010d2,21.0.0-13-g1721dae+760e7a6536,21.0.0-13-g3a573fe+768d78a30a,21.0.0-15-g5a7caf0+f21cbc5713,21.0.0-16-g0fb55c1+b60e2d390c,21.0.0-19-g4cded4ca+71a93a33c0,21.0.0-2-g103fe59+bb20972958,21.0.0-2-g45278ab+04719a4bac,21.0.0-2-g5242d73+3ad5d60fb1,21.0.0-2-g7f82c8f+8babb168e8,21.0.0-2-g8f08a60+06509c8b61,21.0.0-2-g8faa9b5+616205b9df,21.0.0-2-ga326454+8babb168e8,21.0.0-2-gde069b7+5e4aea9c2f,21.0.0-2-gecfae73+1d3a86e577,21.0.0-2-gfc62afb+3ad5d60fb1,21.0.0-25-g1d57be3cd+e73869a214,21.0.0-3-g357aad2+ed88757d29,21.0.0-3-g4a4ce7f+3ad5d60fb1,21.0.0-3-g4be5c26+3ad5d60fb1,21.0.0-3-g65f322c+e0b24896a3,21.0.0-3-g7d9da8d+616205b9df,21.0.0-3-ge02ed75+a9c1acf22d,21.0.0-4-g591bb35+a9c1acf22d,21.0.0-4-g65b4814+b60e2d390c,21.0.0-4-gccdca77+0de219a2bc,21.0.0-4-ge8a399c+6c55c39e83,21.0.0-5-gd00fb1e+05fce91b99,21.0.0-6-gc675373+3ad5d60fb1,21.0.0-64-g1122c245+4fb2b8f86e,21.0.0-7-g04766d7+cd19d05db2,21.0.0-7-gdf92d54+04719a4bac,21.0.0-8-g5674e7b+d1bd76f71f,master-gac4afde19b+a9c1acf22d,w.2021.13
LSST Data Management Base Package
shapeletFunction.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 
24 #include "ndarray/pybind11.h"
25 
27 
28 namespace py = pybind11;
29 using namespace pybind11::literals;
30 
31 namespace lsst {
32 namespace shapelet {
33 
34 PYBIND11_MODULE(shapeletFunction, mod) {
35  py::module::import("lsst.geom");
36  py::module::import("lsst.afw.geom");
37  py::module::import("lsst.afw.image");
38 
39  py::class_<ShapeletFunction, std::shared_ptr<ShapeletFunction>> clsShapeletFunction(mod,
40  "ShapeletFunction");
41 
42  clsShapeletFunction.def_readonly_static("FLUX_FACTOR", &ShapeletFunction::FLUX_FACTOR);
43 
44  clsShapeletFunction.def(py::init<int, BasisTypeEnum>(), "order"_a, "basisType"_a);
45  clsShapeletFunction.def(py::init<int, BasisTypeEnum, ndarray::Array<double, 1, 1> const &>(), "order"_a,
46  "basisType"_a, "coefficients"_a);
47  clsShapeletFunction.def(py::init<int, BasisTypeEnum, double, geom::Point2D const &>(), "order"_a,
48  "basisType"_a, "radius"_a, "center"_a = geom::Point2D());
49  clsShapeletFunction.def(py::init<int, BasisTypeEnum, double, geom::Point2D const &,
50  ndarray::Array<double, 1, 1> const &>(),
51  "order"_a, "basisType"_a, "radius"_a, "center"_a, "coefficients"_a);
52  clsShapeletFunction.def(py::init<int, BasisTypeEnum, afw::geom::ellipses::Ellipse const &>(), "order"_a,
53  "basisType"_a, "ellipse"_a);
54  clsShapeletFunction.def(py::init<int, BasisTypeEnum, afw::geom::ellipses::Ellipse const &,
55  ndarray::Array<double const, 1, 1> const &>(),
56  "order"_a, "basisType"_a, "ellipse"_a, "coefficients"_a);
57  clsShapeletFunction.def(py::init<ShapeletFunction>());
58 
59  clsShapeletFunction.def("getOrder", &ShapeletFunction::getOrder);
60  clsShapeletFunction.def("getEllipse", (afw::geom::ellipses::Ellipse & (ShapeletFunction::*)()) &
61  ShapeletFunction::getEllipse,
62  py::return_value_policy::reference_internal);
63  clsShapeletFunction.def("setEllipse", &ShapeletFunction::setEllipse);
64  clsShapeletFunction.def("getBasisType", &ShapeletFunction::getBasisType);
65  clsShapeletFunction.def("changeBasisType", &ShapeletFunction::changeBasisType);
66  clsShapeletFunction.def("normalize", &ShapeletFunction::normalize, "value"_a = 1.0);
67  clsShapeletFunction.def("getCoefficients", (ndarray::Array<double, 1, 1> const (ShapeletFunction::*)()) &
68  ShapeletFunction::getCoefficients);
69  clsShapeletFunction.def("convolve", &ShapeletFunction::convolve);
70  clsShapeletFunction.def("evaluate", &ShapeletFunction::evaluate);
71  clsShapeletFunction.def("shiftInPlace", &ShapeletFunction::shiftInPlace);
72  clsShapeletFunction.def("transformInPlace", &ShapeletFunction::transformInPlace);
73 
74  py::class_<ShapeletFunctionEvaluator, std::shared_ptr<ShapeletFunctionEvaluator>>
75  clsShapeletFunctionEvaluator(mod, "ShapeletFunctionEvaluator");
76 
77  clsShapeletFunctionEvaluator.def(py::init<ShapeletFunction const &>(), "function"_a);
78 
79  clsShapeletFunctionEvaluator.def("__call__",
80  (double (ShapeletFunctionEvaluator::*)(double, double) const) &
81  ShapeletFunctionEvaluator::operator());
82  clsShapeletFunctionEvaluator.def(
83  "__call__", (double (ShapeletFunctionEvaluator::*)(geom::Point2D const &) const) &
84  ShapeletFunctionEvaluator::operator());
85  clsShapeletFunctionEvaluator.def(
86  "__call__", (double (ShapeletFunctionEvaluator::*)(geom::Extent2D const &) const) &
87  ShapeletFunctionEvaluator::operator());
88  clsShapeletFunctionEvaluator.def("__call__", (ndarray::Array<double, 1, 1> (ShapeletFunctionEvaluator::*)(
89  ndarray::Array<double const, 1> const &,
90  ndarray::Array<double const, 1> const &) const) &
91  ShapeletFunctionEvaluator::operator());
92 
93  clsShapeletFunctionEvaluator.def(
94  "addToImage", (void (ShapeletFunctionEvaluator::*)(ndarray::Array<double, 2, 1> const &,
95  geom::Point2I const &) const) &
97  "array"_a, "xy0"_a = geom::Point2D());
98  clsShapeletFunctionEvaluator.def(
99  "addToImage", (void (ShapeletFunctionEvaluator::*)(afw::image::Image<double> &) const) &
101  "image"_a);
102  clsShapeletFunctionEvaluator.def("integrate", &ShapeletFunctionEvaluator::integrate);
103  clsShapeletFunctionEvaluator.def("computeMoments", &ShapeletFunctionEvaluator::computeMoments);
104  clsShapeletFunctionEvaluator.def("update", &ShapeletFunctionEvaluator::update);
105 }
106 
107 } // shapelet
108 } // lsst
An ellipse defined by an arbitrary BaseCore and a center point.
Definition: Ellipse.h:51
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:58
A 2-d function defined by an expansion onto a Gauss-Laguerre or Gauss-Hermite basis.
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
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.
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
PYBIND11_MODULE(shapeletFunction, mod)
def init()
Definition: tests.py:59
A base class for image defects.