23 #include <pybind11/pybind11.h>
30 using namespace py::literals;
37 template <
typename OutImageT,
typename InImageT,
typename KernelT>
39 mod.def(
"convolve", (
void (*)(OutImageT &, InImageT
const &, KernelT
const &,
40 ConvolutionControl
const &))convolve<OutImageT, InImageT, KernelT>,
41 "convolvedImage"_a,
"inImage"_a,
"kernel"_a,
"convolutionControl"_a = ConvolutionControl());
42 mod.def(
"convolve", (
void (*)(OutImageT &, InImageT
const &, KernelT
const &,
bool,
43 bool))convolve<OutImageT, InImageT, KernelT>,
44 "convolvedImage"_a,
"inImage"_a,
"kernel"_a,
"doNormalize"_a,
"doCopyEdge"_a =
false);
47 template <
typename ImageType1,
typename ImageType2>
50 (
void (*)(ImageType1 &,
double, ImageType2
const &,
double, ImageType2
const &))
scaledPlus);
53 template <
typename ImageType1,
typename ImageType2>
55 declareConvolve<ImageType1, ImageType2, AnalyticKernel>(mod);
56 declareConvolve<ImageType1, ImageType2, DeltaFunctionKernel>(mod);
57 declareConvolve<ImageType1, ImageType2, FixedKernel>(mod);
58 declareConvolve<ImageType1, ImageType2, LinearCombinationKernel>(mod);
59 declareConvolve<ImageType1, ImageType2, SeparableKernel>(mod);
60 declareConvolve<ImageType1, ImageType2, Kernel>(mod);
61 declareScaledPlus<ImageType1, ImageType2>(mod);
64 template <
typename PixelType1,
typename PixelType2>
70 declareByType<M1, M2>(mod);
75 py::class_<ConvolutionControl, std::shared_ptr<ConvolutionControl>> clsConvolutionControl(
76 mod,
"ConvolutionControl");
78 clsConvolutionControl.def(py::init<bool, bool, int>(),
"doNormalize"_a =
true,
"doCopyEdge"_a =
false,
79 "maxInterpolationDistance"_a = 10);
81 clsConvolutionControl.def(
"getDoNormalize", &ConvolutionControl::getDoNormalize);
82 clsConvolutionControl.def(
"getDoCopyEdge", &ConvolutionControl::getDoCopyEdge);
83 clsConvolutionControl.def(
"getMaxInterpolationDistance",
84 &ConvolutionControl::getMaxInterpolationDistance);
85 clsConvolutionControl.def(
"setDoNormalize", &ConvolutionControl::setDoNormalize);
86 clsConvolutionControl.def(
"setDoCopyEdge", &ConvolutionControl::setDoCopyEdge);
87 clsConvolutionControl.def(
"setMaxInterpolationDistance",
88 &ConvolutionControl::setMaxInterpolationDistance);
90 declareAll<double, double>(mod);
91 declareAll<double, float>(mod);
92 declareAll<double, int>(mod);
93 declareAll<double, std::uint16_t>(mod);
94 declareAll<float, float>(mod);
95 declareAll<float, int>(mod);
96 declareAll<float, std::uint16_t>(mod);
97 declareAll<int, int>(mod);
98 declareAll<std::uint16_t, std::uint16_t>(mod);
A class to represent a 2-dimensional array of pixels.
A class to manipulate images, masks, and variance as a single object.
PYBIND11_MODULE(convolveImage, mod)
void scaledPlus(OutImageT &outImage, double c1, InImageT const &inImage1, double c2, InImageT const &inImage2)
Compute the scaled sum of two images.
A base class for image defects.