Loading [MathJax]/extensions/tex2jax.js
LSST Applications g04a91732dc+a777afbe81,g07dc498a13+7e3c5f68a2,g12483e3c20+0145ec33cd,g1409bbee79+7e3c5f68a2,g1a7e361dbc+7e3c5f68a2,g1fd858c14a+9f35e23ec3,g35bb328faa+fcb1d3bbc8,g3ad4f90e5c+0145ec33cd,g3bd4b5ce2c+cbf1bea503,g4e0f332c67+5d362be553,g53246c7159+fcb1d3bbc8,g5477a8d5ce+db04660fe6,g60b5630c4e+0145ec33cd,g623d845a50+0145ec33cd,g6f0c2978f1+3526b51a37,g75b6c65c88+d54b601591,g78460c75b0+2f9a1b4bcd,g786e29fd12+cf7ec2a62a,g7b71ed6315+fcb1d3bbc8,g8852436030+4639f750a5,g89139ef638+7e3c5f68a2,g9125e01d80+fcb1d3bbc8,g919ac25b3e+6220c5324a,g95236ca021+f7a31438ed,g989de1cb63+7e3c5f68a2,g9f33ca652e+2d6fa11d35,gaaedd4e678+7e3c5f68a2,gabe3b4be73+1e0a283bba,gb1101e3267+4a428ef779,gb4a253aaf5+0122250889,gb58c049af0+f03b321e39,gc99c83e5f0+76d20ab76d,gcf25f946ba+4639f750a5,gd6cbbdb0b4+c8606af20c,gde0f65d7ad+3d8a3b7e46,ge278dab8ac+932305ba37,gf795337580+03b96afe58,gfba249425e+fcb1d3bbc8,w.2025.08
LSST Data Management Base Package
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
detect_pybind11.cc File Reference
#include <pybind11/pybind11.h>
#include <pybind11/numpy.h>
#include <pybind11/stl.h>
#include <pybind11/eigen.h>
#include <math.h>
#include <algorithm>
#include <stack>
#include <queue>
#include <vector>
#include <utility>
#include <stdexcept>
#include <iostream>

Go to the source code of this file.

Classes

class  Peak
 A Peak in a Footprint This class is meant to keep track of both the position and flux at the location of a maximum in a Footprint. More...
 
class  Footprint
 

Typedefs

typedef Eigen::Array< int, 4, 1 > Bounds
 
typedef Eigen::Matrix< bool, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > MatrixB
 

Functions

template<typename M>
void get_connected_pixels (const int start_i, const int start_j, py::EigenDRef< const M > image, py::EigenDRef< MatrixB > unchecked, py::EigenDRef< MatrixB > footprint, Eigen::Ref< Bounds > bounds, const double thresh=0)
 
template<typename M>
MatrixB get_connected_multipeak (py::EigenDRef< const M > image, const std::vector< std::vector< int > > &centers, const double thresh=0)
 Proximal operator to trim pixels not connected to one of the source centers.
 
bool sortBrightness (Peak a, Peak b)
 Sort two peaks, placing the brightest peak first.
 
template<typename M>
std::vector< Peakget_peaks (M &image, const double min_separation, const double peak_thresh, const int y0, const int x0)
 
template<typename M>
void maskImage (py::EigenDRef< M > image, py::EigenDRef< MatrixB > footprint)
 
template<typename M, typename P>
std::vector< Footprintget_footprints (py::EigenDRef< const M > image, const double min_separation, const int min_area, const double peak_thresh, const double footprint_thresh, const bool find_peaks=true, const int y0=0, const int x0=0)
 Get all footprints in an image.
 
 PYBIND11_MODULE (detect_pybind11, mod)
 

Typedef Documentation

◆ Bounds

typedef Eigen::Array<int, 4, 1> Bounds

Definition at line 17 of file detect_pybind11.cc.

◆ MatrixB

typedef Eigen::Matrix<bool, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> MatrixB

Definition at line 18 of file detect_pybind11.cc.

Function Documentation

◆ get_connected_multipeak()

template<typename M>
MatrixB get_connected_multipeak ( py::EigenDRef< const M > image,
const std::vector< std::vector< int > > & centers,
const double thresh = 0 )

Proximal operator to trim pixels not connected to one of the source centers.

Definition at line 79 of file detect_pybind11.cc.

83 {
84 const int height = image.rows();
85 const int width = image.cols();
86 MatrixB footprint = MatrixB::Zero(height, width);
88
89 // Seed the queue with peaks
90 for(const auto& center : centers){
91 const int y = center[0];
92 const int x = center[1];
93 if (!footprint(y, x) && image(y, x) > thresh) {
94 footprint(y, x) = true;
95 pixel_queue.emplace(y, x);
96 }
97 }
98
99 // 4-connectivity offsets
100 const std::vector<std::pair<int, int>> offsets = {{-1, 0}, {1, 0}, {0, -1}, {0, 1}};
101
102 // Flood fill
103 while (!pixel_queue.empty()) {
104 auto [i, j] = pixel_queue.front();
105 pixel_queue.pop();
106
107 for (const auto& [di, dj] : offsets) {
108 int ni = i + di;
109 int nj = j + dj;
110 if (ni >= 0 && ni < height && nj >= 0 && nj < width) {
111 if (!footprint(ni, nj) && image(ni, nj) > thresh) {
112 footprint(ni, nj) = true;
113 pixel_queue.emplace(ni, nj);
114 }
115 }
116 }
117 }
118
119 return footprint;
120}
Eigen::Matrix< bool, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > MatrixB
T emplace(T... args)
T empty(T... args)
T front(T... args)
T pop(T... args)

◆ get_connected_pixels()

template<typename M>
void get_connected_pixels ( const int start_i,
const int start_j,
py::EigenDRef< const M > image,
py::EigenDRef< MatrixB > unchecked,
py::EigenDRef< MatrixB > footprint,
Eigen::Ref< Bounds > bounds,
const double thresh = 0 )

Definition at line 24 of file detect_pybind11.cc.

32 {
34 stack.push(std::make_pair(start_i, start_j));
35
36 while (!stack.empty()) {
37 int i, j;
38 std::tie(i, j) = stack.top();
39 stack.pop();
40
41 if (!unchecked(i, j)) {
42 continue;
43 }
44 unchecked(i, j) = false;
45
46 if (image(i, j) > thresh) {
47 footprint(i, j) = true;
48
49 if (i < bounds[0]) {
50 bounds[0] = i;
51 } else if (i > bounds[1]) {
52 bounds[1] = i;
53 }
54 if (j < bounds[2]) {
55 bounds[2] = j;
56 } else if (j > bounds[3]) {
57 bounds[3] = j;
58 }
59
60 if (i > 0 && unchecked(i-1, j)) {
61 stack.push(std::make_pair(i-1, j));
62 }
63 if (i < image.rows() - 1 && unchecked(i+1, j)) {
64 stack.push(std::make_pair(i+1, j));
65 }
66 if (j > 0 && unchecked(i, j-1)) {
67 stack.push(std::make_pair(i, j-1));
68 }
69 if (j < image.cols() - 1 && unchecked(i, j+1)) {
70 stack.push(std::make_pair(i, j+1));
71 }
72 }
73 }
74}
T make_pair(T... args)
T push(T... args)
T tie(T... args)
T top(T... args)

◆ get_footprints()

template<typename M, typename P>
std::vector< Footprint > get_footprints ( py::EigenDRef< const M > image,
const double min_separation,
const int min_area,
const double peak_thresh,
const double footprint_thresh,
const bool find_peaks = true,
const int y0 = 0,
const int x0 = 0 )

Get all footprints in an image.

Parameters
imageThe image to search for footprints
min_separationThe minimum separation (in pixels) between peaks in a footprint
min_areaThe minimum area of a footprint in pixels
peak_threshThe minimum flux of a peak to be detected.
footprint_threshThe minimum flux of a pixel to be included in a footprint
find_peaksIf True, find peaks in each footprint
y0The y-coordinate of the top-left corner of the image
x0The x-coordinate of the top-left corner of the image
Returns
: A list of Footprints

Definition at line 297 of file detect_pybind11.cc.

306 {
307 const int height = image.rows();
308 const int width = image.cols();
309
311 MatrixB unchecked = MatrixB::Ones(height, width);
312 MatrixB footprint = MatrixB::Zero(height, width);
313
314 for(int i=0; i<height; i++){
315 for(int j=0; j<width; j++){
316 Bounds bounds; bounds << i, i, j, j;
317 get_connected_pixels(i, j, image, unchecked, footprint, bounds, footprint_thresh);
318 int subHeight = bounds[1]-bounds[0]+1;
319 int subWidth = bounds[3]-bounds[2]+1;
320 if(subHeight * subWidth > min_area){
321 MatrixB subFootprint = footprint.block(bounds[0], bounds[2], subHeight, subWidth);
322 int area = subFootprint.count();
323 if(area >= min_area){
324 M patch = image.block(bounds[0], bounds[2], subHeight, subWidth);
325 maskImage<M>(patch, subFootprint);
326 std::vector<Peak> _peaks;
327 if(find_peaks){
328 _peaks = get_peaks(
329 patch,
330 min_separation,
331 peak_thresh,
332 bounds[0] + y0,
333 bounds[2] + x0
334 );
335 }
336 // Only add footprints that have at least one peak above the
337 // minimum peak_thresh.
338 if(!_peaks.empty() || !find_peaks){
339 Bounds trueBounds; trueBounds << bounds[0] + y0,
340 bounds[1] + y0, bounds[2] + x0, bounds[3] + x0;
341 footprints.push_back(Footprint(subFootprint, _peaks, trueBounds));
342 }
343 }
344 }
345 footprint.block(bounds[0], bounds[2], subHeight, subWidth) = MatrixB::Zero(subHeight, subWidth);
346 }
347 }
348 return footprints;
349}
Eigen::Array< int, 4, 1 > Bounds
void get_connected_pixels(const int start_i, const int start_j, py::EigenDRef< const M > image, py::EigenDRef< MatrixB > unchecked, py::EigenDRef< MatrixB > footprint, Eigen::Ref< Bounds > bounds, const double thresh=0)
std::vector< Peak > get_peaks(M &image, const double min_separation, const double peak_thresh, const int y0, const int x0)
void maskImage(py::EigenDRef< M > image, py::EigenDRef< MatrixB > footprint)

◆ get_peaks()

template<typename M>
std::vector< Peak > get_peaks ( M & image,
const double min_separation,
const double peak_thresh,
const int y0,
const int x0 )

Sort the peaks in the footprint so that the brightest are first

Definition at line 166 of file detect_pybind11.cc.

172 {
173 const int height = image.rows();
174 const int width = image.cols();
175
176 std::vector<Peak> peaks;
177
178 for(int i=0; i<height; i++){
179 for(int j=0; j<width; j++){
180 if(image(i, j) < peak_thresh){
181 continue;
182 }
183 if(i > 0 && image(i, j) <= image(i-1, j)){
184 continue;
185 }
186 if(i < height-1 && image(i,j) <= image(i+1, j)){
187 continue;
188 }
189 if(j > 0 && image(i, j) <= image(i, j-1)){
190 continue;
191 }
192 if(j < width-1 && image(i,j) <= image(i, j+1)){
193 continue;
194 }
195
196 if(i > 0 && j > 0 && image(i, j) <= image(i-1, j-1)){
197 continue;
198 }
199 if(i < height-1 && j < width-1 && image(i,j) <= image(i+1, j+1)){
200 continue;
201 }
202 if(i < height-1 && j > 0 && image(i, j) <= image(i+1, j-1)){
203 continue;
204 }
205 if(i > 0 && j < width-1 && image(i,j) <= image(i-1, j+1)){
206 continue;
207 }
208
209 peaks.push_back(Peak(i+y0, j+x0, static_cast<double>(image(i, j))));
210 }
211 }
212
213 if(peaks.empty()){
214 return peaks;
215 }
216
218 std::sort (peaks.begin(), peaks.end(), sortBrightness);
219
220 // Remove peaks within min_separation
221 double min_separation2 = min_separation * min_separation;
222 for (size_t i = 0; i < peaks.size() - 1; ++i) {
223 for (size_t j = i + 1; j < peaks.size();) {
224 Peak *p1 = &peaks[i];
225 Peak *p2 = &peaks[j];
226 double dy = p1->getY() - p2->getY();
227 double dx = p1->getX() - p2->getX();
228 double separation2 = dy*dy + dx*dx;
229 if (separation2 < min_separation2) {
230 peaks.erase(peaks.begin() + j);
231 } else {
232 ++j;
233 }
234 }
235 }
236 return peaks;
237}
T begin(T... args)
A Peak in a Footprint This class is meant to keep track of both the position and flux at the location...
bool sortBrightness(Peak a, Peak b)
Sort two peaks, placing the brightest peak first.
T end(T... args)
T erase(T... args)
T push_back(T... args)
T size(T... args)
T sort(T... args)

◆ maskImage()

template<typename M>
void maskImage ( py::EigenDRef< M > image,
py::EigenDRef< MatrixB > footprint )

Definition at line 266 of file detect_pybind11.cc.

269 {
270 const int height = image.rows();
271 const int width = image.cols();
272
273 for(int i=0; i<height; i++){
274 for(int j=0; j<width; j++){
275 if(!footprint(i,j)){
276 image(i,j) = 0;
277 }
278 }
279 }
280}

◆ PYBIND11_MODULE()

PYBIND11_MODULE ( detect_pybind11 ,
mod  )

Definition at line 353 of file detect_pybind11.cc.

353 {
354 typedef Eigen::Matrix<float, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> MatrixF;
355 typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> MatrixD;
356
357 mod.doc() = "Fast detection algorithms implemented in C++";
358
359 mod.def("get_connected_pixels", &get_connected_pixels<MatrixB>,
360 "Create a boolean mask for pixels that are connected");
361 mod.def("get_connected_pixels", &get_connected_pixels<MatrixF>,
362 "Create a boolean mask for pixels that are connected");
363 mod.def("get_connected_pixels", &get_connected_pixels<MatrixD>,
364 "Create a boolean mask for pixels that are connected");
365
366 mod.def("get_connected_multipeak", &get_connected_multipeak<MatrixB>,
367 "Trim pixels not conencted to a center from a list of centers");
368 mod.def("get_connected_multipeak", &get_connected_multipeak<MatrixF>,
369 "Trim pixels not conencted to a center from a list of centers");
370 mod.def("get_connected_multipeak", &get_connected_multipeak<MatrixD>,
371 "Trim pixels not conencted to a center from a list of centers");
372
373 mod.def("get_peaks", &get_peaks<MatrixF>,
374 "Get a list of peaks in a footprint created by get_connected_pixels");
375 mod.def("get_peaks", &get_peaks<MatrixD>,
376 "Get a list of peaks in a footprint created by get_connected_pixels");
377
378 mod.def("get_footprints", &get_footprints<MatrixF, float>,
379 "Create a list of all of the footprints in an image, with their peaks"
380 "image"_a, "min_separation"_a, "min_area"_a, "peak_thresh"_a, "footprint_thresh"_a,
381 "find_peaks"_a=true, "y0"_a=0, "x0"_a=0);
382 mod.def("get_footprints", &get_footprints<MatrixD, double>,
383 "Create a list of all of the footprints in an image, with their peaks"
384 "image"_a, "min_separation"_a, "min_area"_a, "peak_thresh"_a, "footprint_thresh"_a,
385 "find_peaks"_a=true, "y0"_a=0, "x0"_a=0);
386
387 py::class_<Footprint>(mod, "Footprint")
388 .def(py::init<MatrixB, std::vector<Peak>, Bounds>(),
389 "footprint"_a, "peaks"_a, "bounds"_a)
390 .def_property_readonly("data", &Footprint::getFootprint)
391 .def_readwrite("peaks", &Footprint::peaks)
392 .def_property_readonly("bounds", &Footprint::getBounds);
393
394 py::class_<Peak>(mod, "Peak")
395 .def(py::init<int, int, double>(),
396 "y"_a, "x"_a, "flux"_a)
397 .def_property_readonly("y", &Peak::getY)
398 .def_property_readonly("x", &Peak::getX)
399 .def_property_readonly("flux", &Peak::getFlux);
400}
std::vector< Peak > peaks
Bounds getBounds()
MatrixB getFootprint()
double getFlux()
MatrixB get_connected_multipeak(py::EigenDRef< const M > image, const std::vector< std::vector< int > > &centers, const double thresh=0)
Proximal operator to trim pixels not connected to one of the source centers.
std::vector< Footprint > get_footprints(py::EigenDRef< const M > image, const double min_separation, const int min_area, const double peak_thresh, const double footprint_thresh, const bool find_peaks=true, const int y0=0, const int x0=0)
Get all footprints in an image.

◆ sortBrightness()

bool sortBrightness ( Peak a,
Peak b )

Sort two peaks, placing the brightest peak first.

Definition at line 155 of file detect_pybind11.cc.

155 {
156 return a.getFlux() > b.getFlux();
157}