LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
BasicConvolve.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
25 /*
26  * Definition of basicConvolve and convolveWithBruteForce functions declared in detail/ConvolveImage.h
27  */
28 #include <algorithm>
29 #include <cmath>
30 #include <cstdint>
31 #include <sstream>
32 #include <vector>
33 
34 #include "lsst/pex/exceptions.h"
35 #include "lsst/log/Log.h"
36 #include "lsst/geom.h"
39 #include "lsst/afw/math/Kernel.h"
41 
43 
44 namespace {
45 
54 template <typename OutImageT, typename InImageT>
55 void assertDimensionsOK(OutImageT const& convolvedImage, InImageT const& inImage,
56  lsst::afw::math::Kernel const& kernel) {
57  if (convolvedImage.getDimensions() != inImage.getDimensions()) {
59  os << "convolvedImage dimensions = ( " << convolvedImage.getWidth() << ", "
60  << convolvedImage.getHeight() << ") != (" << inImage.getWidth() << ", " << inImage.getHeight()
61  << ") = inImage dimensions";
63  }
64  if (inImage.getWidth() < kernel.getWidth() || inImage.getHeight() < kernel.getHeight()) {
66  os << "inImage dimensions = ( " << inImage.getWidth() << ", " << inImage.getHeight()
67  << ") smaller than (" << kernel.getWidth() << ", " << kernel.getHeight()
68  << ") = kernel dimensions in width and/or height";
70  }
71  if ((kernel.getWidth() < 1) || (kernel.getHeight() < 1)) {
73  os << "kernel dimensions = ( " << kernel.getWidth() << ", " << kernel.getHeight()
74  << ") smaller than (1, 1) in width and/or height";
76  }
77 }
119 template <typename OutPixelT, typename ImageIterT, typename KernelIterT, typename KernelPixelT>
120 inline OutPixelT kernelDotProduct(
121  ImageIterT imageIter,
122  KernelIterT kernelIter,
123  int kWidth)
124 {
125  OutPixelT outPixel(0);
126  for (int x = 0; x < kWidth; ++x, ++imageIter, ++kernelIter) {
127  KernelPixelT kVal = *kernelIter;
128  if (kVal != 0) {
129  outPixel += static_cast<OutPixelT>((*imageIter) * kVal);
130  }
131  }
132  return outPixel;
133 }
134 } // anonymous namespace
135 
136 namespace lsst {
137 namespace afw {
138 namespace math {
139 namespace detail {
140 
141 template <typename OutImageT, typename InImageT>
142 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage, math::Kernel const& kernel,
143  math::ConvolutionControl const& convolutionControl) {
144  // Because convolve isn't a method of Kernel we can't always use Kernel's vtbl to dynamically
145  // dispatch the correct version of basicConvolve. The case that fails is convolving with a kernel
146  // obtained from a pointer or reference to a Kernel (base class), e.g. as used in LinearCombinationKernel.
147  if (IS_INSTANCE(kernel, math::DeltaFunctionKernel)) {
148  LOGL_DEBUG("TRACE3.afw.math.convolve.basicConvolve",
149  "generic basicConvolve: dispatch to DeltaFunctionKernel basicConvolve");
150  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::DeltaFunctionKernel const*>(&kernel),
151  convolutionControl);
152  return;
153  } else if (IS_INSTANCE(kernel, math::SeparableKernel)) {
154  LOGL_DEBUG("TRACE3.afw.math.convolve.basicConvolve",
155  "generic basicConvolve: dispatch to SeparableKernel basicConvolve");
156  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::SeparableKernel const*>(&kernel),
157  convolutionControl);
158  return;
159  } else if (IS_INSTANCE(kernel, math::LinearCombinationKernel) && kernel.isSpatiallyVarying()) {
160  LOGL_DEBUG(
161  "TRACE3.afw.math.convolve.basicConvolve",
162  "generic basicConvolve: dispatch to spatially varying LinearCombinationKernel basicConvolve");
163  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::LinearCombinationKernel const*>(&kernel),
164  convolutionControl);
165  return;
166  }
167  // OK, use general (and slower) form
168  if (kernel.isSpatiallyVarying() && (convolutionControl.getMaxInterpolationDistance() > 1)) {
169  // use linear interpolation
170  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
171  "generic basicConvolve: using linear interpolation");
172  convolveWithInterpolation(convolvedImage, inImage, kernel, convolutionControl);
173 
174  } else {
175  // use brute force
176  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve", "generic basicConvolve: using brute force");
177  convolveWithBruteForce(convolvedImage, inImage, kernel, convolutionControl);
178  }
179 }
180 
181 template <typename OutImageT, typename InImageT>
182 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage,
183  math::DeltaFunctionKernel const& kernel,
184  math::ConvolutionControl const& convolutionControl) {
185  assert(!kernel.isSpatiallyVarying());
186  assertDimensionsOK(convolvedImage, inImage, kernel);
187 
188  int const mImageWidth = inImage.getWidth(); // size of input region
189  int const mImageHeight = inImage.getHeight();
190  int const cnvWidth = mImageWidth + 1 - kernel.getWidth();
191  int const cnvHeight = mImageHeight + 1 - kernel.getHeight();
192  int const cnvStartX = kernel.getCtrX();
193  int const cnvStartY = kernel.getCtrY();
194  int const inStartX = kernel.getPixel().getX();
195  int const inStartY = kernel.getPixel().getY();
196 
197  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve", "DeltaFunctionKernel basicConvolve");
198 
199  for (int i = 0; i < cnvHeight; ++i) {
200  typename InImageT::x_iterator inPtr = inImage.x_at(inStartX, i + inStartY);
201  for (typename OutImageT::x_iterator cnvPtr = convolvedImage.x_at(cnvStartX, i + cnvStartY),
202  cnvEnd = cnvPtr + cnvWidth;
203  cnvPtr != cnvEnd; ++cnvPtr, ++inPtr) {
204  *cnvPtr = *inPtr;
205  }
206  }
207 }
208 
209 template <typename OutImageT, typename InImageT>
210 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage,
211  math::LinearCombinationKernel const& kernel,
212  math::ConvolutionControl const& convolutionControl) {
213  if (!kernel.isSpatiallyVarying()) {
214  // use the standard algorithm for the spatially invariant case
215  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
216  "basicConvolve for LinearCombinationKernel: spatially invariant; using brute force");
217  return convolveWithBruteForce(convolvedImage, inImage, kernel, convolutionControl.getDoNormalize());
218  } else {
219  // refactor the kernel if this is reasonable and possible;
220  // then use the standard algorithm for the spatially varying case
221  std::shared_ptr<Kernel> refKernelPtr; // possibly refactored version of kernel
222  if (static_cast<int>(kernel.getNKernelParameters()) > kernel.getNSpatialParameters()) {
223  // refactoring will speed convolution, so try it
224  refKernelPtr = kernel.refactor();
225  if (!refKernelPtr) {
226  refKernelPtr = kernel.clone();
227  }
228  } else {
229  // too few basis kernels for refactoring to be worthwhile
230  refKernelPtr = kernel.clone();
231  }
232  if (convolutionControl.getMaxInterpolationDistance() > 1) {
233  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
234  "basicConvolve for LinearCombinationKernel: using interpolation");
235  return convolveWithInterpolation(convolvedImage, inImage, *refKernelPtr, convolutionControl);
236  } else {
237  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
238  "basicConvolve for LinearCombinationKernel: maxInterpolationError < 0; using brute "
239  "force");
240  return convolveWithBruteForce(convolvedImage, inImage, *refKernelPtr,
241  convolutionControl.getDoNormalize());
242  }
243  }
244 }
245 
246 template <typename OutImageT, typename InImageT>
247 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage, math::SeparableKernel const& kernel,
248  math::ConvolutionControl const& convolutionControl) {
249  typedef typename math::Kernel::Pixel KernelPixel;
250  typedef typename std::vector<KernelPixel> KernelVector;
251  typedef KernelVector::const_iterator KernelIterator;
252  typedef typename InImageT::const_x_iterator InXIterator;
253  typedef typename InImageT::const_xy_locator InXYLocator;
254  typedef typename OutImageT::x_iterator OutXIterator;
255  typedef typename OutImageT::y_iterator OutYIterator;
256  typedef typename OutImageT::SinglePixel OutPixel;
257 
258  assertDimensionsOK(convolvedImage, inImage, kernel);
259 
260  lsst::geom::Box2I const fullBBox = inImage.getBBox(image::LOCAL);
261  lsst::geom::Box2I const goodBBox = kernel.shrinkBBox(fullBBox);
262 
263  KernelVector kernelXVec(kernel.getWidth());
264  KernelVector kernelYVec(kernel.getHeight());
265 
266  if (kernel.isSpatiallyVarying()) {
267  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
268  "SeparableKernel basicConvolve: kernel is spatially varying");
269 
270  for (int cnvY = goodBBox.getMinY(); cnvY <= goodBBox.getMaxY(); ++cnvY) {
271  double const rowPos = inImage.indexToPosition(cnvY, image::Y);
272 
273  InXYLocator inImLoc = inImage.xy_at(0, cnvY - goodBBox.getMinY());
274  OutXIterator cnvXIter = convolvedImage.row_begin(cnvY) + goodBBox.getMinX();
275  for (int cnvX = goodBBox.getMinX(); cnvX <= goodBBox.getMaxX();
276  ++cnvX, ++inImLoc.x(), ++cnvXIter) {
277  double const colPos = inImage.indexToPosition(cnvX, image::X);
278 
279  KernelPixel kSum = kernel.computeVectors(kernelXVec, kernelYVec,
280  convolutionControl.getDoNormalize(), colPos, rowPos);
281 
282  // why does this trigger warnings? It did not in the past.
283  *cnvXIter = math::convolveAtAPoint<OutImageT, InImageT>(inImLoc, kernelXVec, kernelYVec);
284  if (convolutionControl.getDoNormalize()) {
285  *cnvXIter = *cnvXIter / kSum;
286  }
287  }
288  }
289  } else {
290  // kernel is spatially invariant
291  // The basic sequence:
292  // - For each output row:
293  // - Compute x-convolved data: a kernel height's strip of input image convolved with kernel x vector
294  // - Compute one row of output by dotting each column of x-convolved data with the kernel y vector
295  // The x-convolved data is stored in a kernel-height by good-width buffer.
296  // This is circular buffer along y (to avoid shifting pixels before setting each new row);
297  // so for each new row the kernel y vector is rotated to match the order of the x-convolved data.
298 
299  LOGL_DEBUG("TRACE2.afw.math.convolve.basicConvolve",
300  "SeparableKernel basicConvolve: kernel is spatially invariant");
301 
302  kernel.computeVectors(kernelXVec, kernelYVec, convolutionControl.getDoNormalize());
303  KernelIterator const kernelXVecBegin = kernelXVec.begin();
304  KernelIterator const kernelYVecBegin = kernelYVec.begin();
305 
306  // buffer for x-convolved data
307  OutImageT buffer(lsst::geom::Extent2I(goodBBox.getWidth(), kernel.getHeight()));
308 
309  // pre-fill x-convolved data buffer with all but one row of data
310  int yInd = 0; // during initial fill bufY = inImageY
311  int const yPrefillEnd = buffer.getHeight() - 1;
312  for (; yInd < yPrefillEnd; ++yInd) {
313  OutXIterator bufXIter = buffer.x_at(0, yInd);
314  OutXIterator const bufXEnd = buffer.x_at(goodBBox.getWidth(), yInd);
315  InXIterator inXIter = inImage.x_at(0, yInd);
316  for (; bufXIter != bufXEnd; ++bufXIter, ++inXIter) {
317  *bufXIter = kernelDotProduct<OutPixel, InXIterator, KernelIterator, KernelPixel>(
318  inXIter, kernelXVecBegin, kernel.getWidth());
319  }
320  }
321 
322  // compute output pixels using the sequence described above
323  int inY = yPrefillEnd;
324  int bufY = yPrefillEnd;
325  int cnvY = goodBBox.getMinY();
326  while (true) {
327  // fill next buffer row and compute output row
328  InXIterator inXIter = inImage.x_at(0, inY);
329  OutXIterator bufXIter = buffer.x_at(0, bufY);
330  OutXIterator cnvXIter = convolvedImage.x_at(goodBBox.getMinX(), cnvY);
331  for (int bufX = 0; bufX < goodBBox.getWidth(); ++bufX, ++cnvXIter, ++bufXIter, ++inXIter) {
332  // note: bufXIter points to the row of the buffer that is being updated,
333  // whereas bufYIter points to row 0 of the buffer
334  *bufXIter = kernelDotProduct<OutPixel, InXIterator, KernelIterator, KernelPixel>(
335  inXIter, kernelXVecBegin, kernel.getWidth());
336 
337  OutYIterator bufYIter = buffer.y_at(bufX, 0);
338  *cnvXIter = kernelDotProduct<OutPixel, OutYIterator, KernelIterator, KernelPixel>(
339  bufYIter, kernelYVecBegin, kernel.getHeight());
340  }
341 
342  // test for done now, instead of the start of the loop,
343  // to avoid an unnecessary extra rotation of the kernel Y vector
344  if (cnvY >= goodBBox.getMaxY()) break;
345 
346  // update y indices, including bufY, and rotate the kernel y vector to match
347  ++inY;
348  bufY = (bufY + 1) % kernel.getHeight();
349  ++cnvY;
350  std::rotate(kernelYVec.begin(), kernelYVec.end() - 1, kernelYVec.end());
351  }
352  }
353 }
354 
355 template <typename OutImageT, typename InImageT>
356 void convolveWithBruteForce(OutImageT& convolvedImage, InImageT const& inImage, math::Kernel const& kernel,
357  math::ConvolutionControl const& convolutionControl) {
358  bool doNormalize = convolutionControl.getDoNormalize();
359 
360  typedef typename math::Kernel::Pixel KernelPixel;
361  typedef image::Image<KernelPixel> KernelImage;
362 
363  typedef typename KernelImage::const_x_iterator KernelXIterator;
364  typedef typename KernelImage::const_xy_locator KernelXYLocator;
365  typedef typename InImageT::const_x_iterator InXIterator;
366  typedef typename InImageT::const_xy_locator InXYLocator;
367  typedef typename OutImageT::x_iterator OutXIterator;
368  typedef typename OutImageT::SinglePixel OutPixel;
369 
370  assertDimensionsOK(convolvedImage, inImage, kernel);
371 
372  int const inImageWidth = inImage.getWidth();
373  int const inImageHeight = inImage.getHeight();
374  int const kWidth = kernel.getWidth();
375  int const kHeight = kernel.getHeight();
376  int const cnvWidth = inImageWidth + 1 - kernel.getWidth();
377  int const cnvHeight = inImageHeight + 1 - kernel.getHeight();
378  int const cnvStartX = kernel.getCtrX();
379  int const cnvStartY = kernel.getCtrY();
380  int const cnvEndX = cnvStartX + cnvWidth; // end index + 1
381  int const cnvEndY = cnvStartY + cnvHeight; // end index + 1
382 
383  KernelImage kernelImage(kernel.getDimensions());
384  KernelXYLocator const kernelLoc = kernelImage.xy_at(0, 0);
385 
386  if (kernel.isSpatiallyVarying()) {
387  LOGL_DEBUG("TRACE4.afw.math.convolve.convolveWithBruteForce",
388  "convolveWithBruteForce: kernel is spatially varying");
389 
390  for (int cnvY = cnvStartY; cnvY != cnvEndY; ++cnvY) {
391  double const rowPos = inImage.indexToPosition(cnvY, image::Y);
392 
393  InXYLocator inImLoc = inImage.xy_at(0, cnvY - cnvStartY);
394  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
395  for (int cnvX = cnvStartX; cnvX != cnvEndX; ++cnvX, ++inImLoc.x(), ++cnvXIter) {
396  double const colPos = inImage.indexToPosition(cnvX, image::X);
397 
398  KernelPixel kSum = kernel.computeImage(kernelImage, false, colPos, rowPos);
399  *cnvXIter = math::convolveAtAPoint<OutImageT, InImageT>(inImLoc, kernelLoc, kWidth, kHeight);
400  if (doNormalize) {
401  *cnvXIter = *cnvXIter / kSum;
402  }
403  }
404  }
405  } else {
406  LOGL_DEBUG("TRACE4.afw.math.convolve.convolveWithBruteForce",
407  "convolveWithBruteForce: kernel is spatially invariant");
408 
409  (void)kernel.computeImage(kernelImage, doNormalize);
410 
411  for (int inStartY = 0, cnvY = cnvStartY; inStartY < cnvHeight; ++inStartY, ++cnvY) {
412  KernelXIterator kernelXIter = kernelImage.x_at(0, 0);
413  InXIterator inXIter = inImage.x_at(0, inStartY);
414  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
415  for (int x = 0; x < cnvWidth; ++x, ++cnvXIter, ++inXIter) {
416  *cnvXIter = kernelDotProduct<OutPixel, InXIterator, KernelXIterator, KernelPixel>(
417  inXIter, kernelXIter, kWidth);
418  }
419  for (int kernelY = 1, inY = inStartY + 1; kernelY < kHeight; ++inY, ++kernelY) {
420  KernelXIterator kernelXIter = kernelImage.x_at(0, kernelY);
421  InXIterator inXIter = inImage.x_at(0, inY);
422  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
423  for (int x = 0; x < cnvWidth; ++x, ++cnvXIter, ++inXIter) {
424  *cnvXIter += kernelDotProduct<OutPixel, InXIterator, KernelXIterator, KernelPixel>(
425  inXIter, kernelXIter, kWidth);
426  }
427  }
428  }
429  }
430 }
431 
432 /*
433  * Explicit instantiation
434  */
436 #define IMAGE(PIXTYPE) image::Image<PIXTYPE>
437 #define MASKEDIMAGE(PIXTYPE) image::MaskedImage<PIXTYPE, image::MaskPixel, image::VariancePixel>
438 #define NL /* */
439 // Instantiate Image or MaskedImage versions
440 #define INSTANTIATE_IM_OR_MI(IMGMACRO, OUTPIXTYPE, INPIXTYPE) \
441  template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, math::Kernel const&, \
442  math::ConvolutionControl const&); \
443  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
444  math::DeltaFunctionKernel const&, math::ConvolutionControl const&); \
445  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
446  math::LinearCombinationKernel const&, math::ConvolutionControl const&); \
447  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
448  math::SeparableKernel const&, math::ConvolutionControl const&); \
449  NL template void convolveWithBruteForce(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
450  math::Kernel const&, math::ConvolutionControl const&);
451 // Instantiate both Image and MaskedImage versions
452 #define INSTANTIATE(OUTPIXTYPE, INPIXTYPE) \
453  INSTANTIATE_IM_OR_MI(IMAGE, OUTPIXTYPE, INPIXTYPE) \
454  INSTANTIATE_IM_OR_MI(MASKEDIMAGE, OUTPIXTYPE, INPIXTYPE)
455 
456 INSTANTIATE(double, double)
457 INSTANTIATE(double, float)
458 INSTANTIATE(double, int)
459 INSTANTIATE(double, std::uint16_t)
460 INSTANTIATE(float, float)
461 INSTANTIATE(float, int)
463 INSTANTIATE(int, int)
466 } // namespace detail
467 } // namespace math
468 } // namespace afw
469 } // namespace lsst
int getHeight() const
Return the Kernel&#39;s height.
Definition: Kernel.h:230
int getCtrX() const
Return x index of kernel&#39;s center.
Definition: Kernel.h:244
int getCtrY() const
Return y index of kernel&#39;s center.
Definition: Kernel.h:255
Parameters to control convolution.
Definition: ConvolveImage.h:50
A kernel described by a pair of functions: func(x, y) = colFunc(x) * rowFunc(y)
Definition: Kernel.h:907
std::shared_ptr< Kernel > clone() const override
Return a pointer to a deep copy of this kernel.
void convolveWithBruteForce(OutImageT &convolvedImage, InImageT const &inImage, lsst::afw::math::Kernel const &kernel, lsst::afw::math::ConvolutionControl const &convolutionControl)
Convolve an Image or MaskedImage with a Kernel by computing the kernel image at every point...
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition: Log.h:504
LSST DM logging module built on log4cxx.
bool isSpatiallyVarying() const
Return true iff the kernel is spatially varying (has a spatial function)
Definition: Kernel.h:380
double computeVectors(std::vector< Pixel > &colList, std::vector< Pixel > &rowList, bool doNormalize, double x=0.0, double y=0.0) const
Compute the column and row arrays in place, where kernel(col, row) = colList(col) * rowList(row) ...
A base class for image defects.
int getMaxY() const noexcept
Definition: Box.h:162
A kernel that is a linear combination of fixed basis kernels.
Definition: Kernel.h:750
T str(T... args)
int getWidth() const noexcept
Definition: Box.h:187
double computeImage(lsst::afw::image::Image< Pixel > &image, bool doNormalize, double x=0.0, double y=0.0) const
Compute an image (pixellized representation of the kernel) in place.
Definition: Kernel.cc:85
void convolveWithInterpolation(OutImageT &outImage, InImageT const &inImage, lsst::afw::math::Kernel const &kernel, ConvolutionControl const &convolutionControl)
Convolve an Image or MaskedImage with a spatially varying Kernel using linear interpolation.
int getMaxX() const noexcept
Definition: Box.h:161
double x
int getWidth() const
Return the Kernel&#39;s width.
Definition: Kernel.h:225
lsst::geom::Extent2I const getDimensions() const
Return the Kernel&#39;s dimensions (width, height)
Definition: Kernel.h:213
unsigned int getNKernelParameters() const
Return the number of kernel parameters (0 if none)
Definition: Kernel.h:269
int getMinX() const noexcept
Definition: Box.h:157
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
STL class.
lsst::geom::Box2I shrinkBBox(lsst::geom::Box2I const &bbox) const
Given a bounding box for an image one wishes to convolve with this kernel, return the bounding box fo...
Definition: Kernel.cc:183
T rotate(T... args)
int getNSpatialParameters() const
Return the number of spatial parameters (0 if not spatially varying)
Definition: Kernel.h:274
Reports invalid arguments.
Definition: Runtime.h:66
lsst::geom::Point2I getPixel() const
Definition: Kernel.h:717
std::ostream * os
Definition: Schema.cc:746
std::shared_ptr< Kernel > refactor() const
Refactor the kernel as a linear combination of N bases where N is the number of parameters for the sp...
Kernels are used for convolution with MaskedImages and (eventually) Images.
Definition: Kernel.h:111
An integer coordinate rectangle.
Definition: Box.h:55
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:58
#define IS_INSTANCE(A, B)
Definition: Convolve.h:40
void basicConvolve(OutImageT &convolvedImage, InImageT const &inImage, lsst::afw::math::Kernel const &kernel, lsst::afw::math::ConvolutionControl const &convolutionControl)
Low-level convolution function that does not set edge pixels.
A kernel that has only one non-zero pixel (of value 1)
Definition: Kernel.h:690
int getMinY() const noexcept
Definition: Box.h:158
#define INSTANTIATE(FROMSYS, TOSYS)
Definition: Detector.cc:484