LSST Applications g063fba187b+fee0456c91,g0f08755f38+ea96e5a5a3,g1653933729+a8ce1bb630,g168dd56ebc+a8ce1bb630,g1a2382251a+90257ff92a,g20f6ffc8e0+ea96e5a5a3,g217e2c1bcf+937a289c59,g28da252d5a+daa7da44eb,g2bbee38e9b+253935c60e,g2bc492864f+253935c60e,g3156d2b45e+6e55a43351,g32e5bea42b+31359a2a7a,g347aa1857d+253935c60e,g35bb328faa+a8ce1bb630,g3a166c0a6a+253935c60e,g3b1af351f3+a8ce1bb630,g3e281a1b8c+c5dd892a6c,g414038480c+416496e02f,g41af890bb2+afe91b1188,g599934f4f4+0db33f7991,g7af13505b9+e36de7bce6,g80478fca09+da231ba887,g82479be7b0+a4516e59e3,g858d7b2824+ea96e5a5a3,g89c8672015+f4add4ffd5,g9125e01d80+a8ce1bb630,ga5288a1d22+bc6ab8dfbd,gb58c049af0+d64f4d3760,gc28159a63d+253935c60e,gcab2d0539d+3f2b72788c,gcf0d15dbbd+4ea9c45075,gda6a2b7d83+4ea9c45075,gdaeeff99f8+1711a396fd,ge79ae78c31+253935c60e,gef2f8181fd+3031e3cf99,gf0baf85859+c1f95f4921,gfa517265be+ea96e5a5a3,gfa999e8aa5+17cd334064,w.2024.50
LSST Data Management Base Package
Loading...
Searching...
No Matches
Static Public Member Functions | List of all members
lsst::meas::algorithms::CloughTocher2DInterpolatorUtils Class Reference

This class contains static utility methods that are used by the CloughTocher2DInterpolatorTask and exists solely to provide a namespace. More...

#include <CloughTocher2DInterpolatorUtils.h>

Static Public Member Functions

static std::pair< ndarray::Array< float, 2, 2 >, ndarray::Array< float, 2, 2 > > findGoodPixelsAroundBadPixels (afw::image::MaskedImage< PixelT, afw::image::MaskPixel, afw::image::VariancePixel > const &mimage, std::vector< std::string > const &maskPlanes, int const buffer)
 Find the positions of bad pixels as defined by the masks, and also find the location and pixel value of good pixels around the bad pixels.
 
static void updateArrayFromImage (ndarray::Array< float, 2, 2 > const &pixelArray, afw::image::Image< PixelT > const &image)
 Update the values (third column) of the pixelArray from the image.
 
static void updateImageFromArray (afw::image::Image< PixelT > &image, ndarray::Array< float const, 2, 2 > const &pixelArray)
 Update the pixel values of the image from the pixelArray.
 

Detailed Description

This class contains static utility methods that are used by the CloughTocher2DInterpolatorTask and exists solely to provide a namespace.

Definition at line 44 of file CloughTocher2DInterpolatorUtils.h.

Member Function Documentation

◆ findGoodPixelsAroundBadPixels()

std::pair< ndarray::Array< float, 2, 2 >, ndarray::Array< float, 2, 2 > > lsst::meas::algorithms::CloughTocher2DInterpolatorUtils::findGoodPixelsAroundBadPixels ( afw::image::MaskedImage< PixelT, afw::image::MaskPixel, afw::image::VariancePixel > const & mimage,
std::vector< std::string > const & maskPlanes,
int const buffer )
static

Find the positions of bad pixels as defined by the masks, and also find the location and pixel value of good pixels around the bad pixels.

Parameters
[in]mimageMaskedImage to find the pixels from.
[in]maskPlanesList of mask planes to consider as bad pixels.
[in]bufferNumber of pixels to dilate the bad pixels by.
Returns
A pair of arrays, the first containing the bad pixels and the second containing the good pixels.
Note
The bad pixels array has shape (N, 3) where N is the number of bad pixels. The first column is the x coordinate, the second column is the y coordinate, and the third column is the pixel value. The good pixels array has shape (M, 3) where M is the number of good pixels and has the same format as the previous one.

Definition at line 43 of file CloughTocher2DInterpolatorUtils.cc.

45 {
46 auto bitMask = afw::image::Mask<>::getPlaneBitMask(maskPlanes);
47
48 auto badSpanSet = afw::geom::SpanSet::fromMask(
49 *mimage.getMask(), [bitMask](afw::image::MaskPixel pixel) { return bool(pixel & bitMask); });
50 ndarray::Array<float, 2, 2> badPixelArray = ndarray::allocate(badSpanSet->getArea(), 3);
51 badSpanSet->applyFunctor(
52 [](geom::Point2I pt, float &x, float &y, float &valueOut, float valueIn) {
53 x = pt.getX();
54 y = pt.getY();
55 valueOut = valueIn;
56 },
57 ndFlat(badPixelArray[ndarray::view()(0)].shallow()),
58 ndFlat(badPixelArray[ndarray::view()(1)].shallow()),
59 ndFlat(badPixelArray[ndarray::view()(2)].shallow()),
60 *mimage.getImage());
61
62 auto allSpanSet = badSpanSet->dilated(buffer, afw::geom::Stencil::BOX);
63 auto goodSpanSet = allSpanSet->intersectNot(*badSpanSet);
64 goodSpanSet = goodSpanSet->clippedTo(mimage.getBBox());
65
66 badSpanSet.reset();
67 allSpanSet.reset();
68
69 ndarray::Array<float, 2, 2> goodPixelArray = ndarray::allocate(goodSpanSet->getArea(), 3);
70 goodSpanSet->applyFunctor(
71 [](geom::Point2I pt, float &x, float &y, float &valueOut, float valueIn) {
72 x = pt.getX();
73 y = pt.getY();
74 valueOut = valueIn;
75 },
76 ndFlat(goodPixelArray[ndarray::view()(0)].shallow()),
77 ndFlat(goodPixelArray[ndarray::view()(1)].shallow()),
78 ndFlat(goodPixelArray[ndarray::view()(2)].shallow()),
79 *mimage.getImage());
80
81 goodSpanSet.reset();
82
83 return std::make_pair(badPixelArray, goodPixelArray);
84}
int y
Definition SpanSet.cc:48
static std::shared_ptr< geom::SpanSet > fromMask(image::Mask< T > const &mask, UnaryPredicate comparator=details::AnyBitSetFunctor< T >())
Create a SpanSet from a mask.
Definition SpanSet.h:644
T make_pair(T... args)
std::int32_t MaskPixel
default type for Masks and MaskedImage Masks
details::FlatNdGetter< T, inA, inB > ndFlat(ndarray::Array< T, inA, inB > const &array)
Marks a ndarray to be interpreted as a 1D vector when applying a functor from a SpanSet.

◆ updateArrayFromImage()

void lsst::meas::algorithms::CloughTocher2DInterpolatorUtils::updateArrayFromImage ( ndarray::Array< float, 2, 2 > const & pixelArray,
afw::image::Image< PixelT > const & image )
static

Update the values (third column) of the pixelArray from the image.

Parameters
[out]pixelArrayThe three-column array to update.
[in]imageThe image to update the pixel values from.
Note
The pixelArray is typically one of the arrays returned by findGoodPixelsAroundBadPixels.
See also
updateImageFromArray

Definition at line 86 of file CloughTocher2DInterpolatorUtils.cc.

87 {
88 afw::image::CheckIndices docheck(true);
89 auto x0 = image.getX0();
90 auto y0 = image.getY0();
91 for (auto row : pixelArray) {
92 int x = row[0] - x0;
93 int y = row[1] - y0;
94 row[2] = image(x, y, docheck);
95 }
96}
afw::table::Key< afw::table::Array< ImagePixelT > > image
int row
Definition CR.cc:145

◆ updateImageFromArray()

void lsst::meas::algorithms::CloughTocher2DInterpolatorUtils::updateImageFromArray ( afw::image::Image< PixelT > & image,
ndarray::Array< float const, 2, 2 > const & pixelArray )
static

Update the pixel values of the image from the pixelArray.

Parameters
[out]imageThe image to update.
[in]pixelArrayThe three-column array to update the pixel values from.
Note
The pixelArray is typically one of the arrays returned by findGoodPixelsAroundBadPixels.
See also
updateArrayFromImage

Definition at line 98 of file CloughTocher2DInterpolatorUtils.cc.

99 {
100 afw::image::CheckIndices docheck(true);
101 auto x0 = image.getX0();
102 auto y0 = image.getY0();
103 for (auto row : pixelArray) {
104 int x = row[0] - x0;
105 int y = row[1] - y0;
106 image(x, y, docheck) = row[2];
107 }
108}

The documentation for this class was generated from the following files: