LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
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 <cstdint>
30 #include <sstream>
31 #include <vector>
32 
33 #include "lsst/pex/exceptions.h"
34 #include "lsst/log/Log.h"
35 #include "lsst/geom.h"
38 #include "lsst/afw/math/Kernel.h"
40 
42 
43 namespace {
44 
53 template <typename OutImageT, typename InImageT>
54 void assertDimensionsOK(OutImageT const& convolvedImage, InImageT const& inImage,
55  lsst::afw::math::Kernel const& kernel) {
56  if (convolvedImage.getDimensions() != inImage.getDimensions()) {
58  os << "convolvedImage dimensions = ( " << convolvedImage.getWidth() << ", "
59  << convolvedImage.getHeight() << ") != (" << inImage.getWidth() << ", " << inImage.getHeight()
60  << ") = inImage dimensions";
62  }
63  if (inImage.getWidth() < kernel.getWidth() || inImage.getHeight() < kernel.getHeight()) {
65  os << "inImage dimensions = ( " << inImage.getWidth() << ", " << inImage.getHeight()
66  << ") smaller than (" << kernel.getWidth() << ", " << kernel.getHeight()
67  << ") = kernel dimensions in width and/or height";
69  }
70  if ((kernel.getWidth() < 1) || (kernel.getHeight() < 1)) {
72  os << "kernel dimensions = ( " << kernel.getWidth() << ", " << kernel.getHeight()
73  << ") smaller than (1, 1) in width and/or height";
75  }
76 }
118 template <typename OutPixelT, typename ImageIterT, typename KernelIterT, typename KernelPixelT>
119 inline OutPixelT kernelDotProduct(
120  ImageIterT imageIter,
121  KernelIterT kernelIter,
122  int kWidth)
123 {
124  OutPixelT outPixel(0);
125  for (int x = 0; x < kWidth; ++x, ++imageIter, ++kernelIter) {
126  KernelPixelT kVal = *kernelIter;
127  if (kVal != 0) {
128  outPixel += static_cast<OutPixelT>((*imageIter) * kVal);
129  }
130  }
131  return outPixel;
132 }
133 } // anonymous namespace
134 
135 namespace lsst {
136 namespace afw {
137 namespace math {
138 namespace detail {
139 
140 template <typename OutImageT, typename InImageT>
141 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage, math::Kernel const& kernel,
142  math::ConvolutionControl const& convolutionControl) {
143  // Because convolve isn't a method of Kernel we can't always use Kernel's vtbl to dynamically
144  // dispatch the correct version of basicConvolve. The case that fails is convolving with a kernel
145  // obtained from a pointer or reference to a Kernel (base class), e.g. as used in LinearCombinationKernel.
146  if (IS_INSTANCE(kernel, math::DeltaFunctionKernel)) {
147  LOGL_DEBUG("TRACE3.lsst.afw.math.convolve.basicConvolve",
148  "generic basicConvolve: dispatch to DeltaFunctionKernel basicConvolve");
149  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::DeltaFunctionKernel const*>(&kernel),
150  convolutionControl);
151  return;
152  } else if (IS_INSTANCE(kernel, math::SeparableKernel)) {
153  LOGL_DEBUG("TRACE3.lsst.afw.math.convolve.basicConvolve",
154  "generic basicConvolve: dispatch to SeparableKernel basicConvolve");
155  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::SeparableKernel const*>(&kernel),
156  convolutionControl);
157  return;
158  } else if (IS_INSTANCE(kernel, math::LinearCombinationKernel) && kernel.isSpatiallyVarying()) {
159  LOGL_DEBUG(
160  "TRACE3.afw.math.convolve.basicConvolve",
161  "generic basicConvolve: dispatch to spatially varying LinearCombinationKernel basicConvolve");
162  basicConvolve(convolvedImage, inImage, *dynamic_cast<math::LinearCombinationKernel const*>(&kernel),
163  convolutionControl);
164  return;
165  }
166  // OK, use general (and slower) form
167  if (kernel.isSpatiallyVarying() && (convolutionControl.getMaxInterpolationDistance() > 1)) {
168  // use linear interpolation
169  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
170  "generic basicConvolve: using linear interpolation");
171  convolveWithInterpolation(convolvedImage, inImage, kernel, convolutionControl);
172 
173  } else {
174  // use brute force
175  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve", "generic basicConvolve: using brute force");
176  convolveWithBruteForce(convolvedImage, inImage, kernel, convolutionControl);
177  }
178 }
179 
180 template <typename OutImageT, typename InImageT>
181 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage,
182  math::DeltaFunctionKernel const& kernel,
183  math::ConvolutionControl const& convolutionControl) {
184  assert(!kernel.isSpatiallyVarying());
185  assertDimensionsOK(convolvedImage, inImage, kernel);
186 
187  int const mImageWidth = inImage.getWidth(); // size of input region
188  int const mImageHeight = inImage.getHeight();
189  int const cnvWidth = mImageWidth + 1 - kernel.getWidth();
190  int const cnvHeight = mImageHeight + 1 - kernel.getHeight();
191  int const cnvStartX = kernel.getCtr().getX();
192  int const cnvStartY = kernel.getCtr().getY();
193  int const inStartX = kernel.getPixel().getX();
194  int const inStartY = kernel.getPixel().getY();
195 
196  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve", "DeltaFunctionKernel basicConvolve");
197 
198  for (int i = 0; i < cnvHeight; ++i) {
199  typename InImageT::x_iterator inPtr = inImage.x_at(inStartX, i + inStartY);
200  for (typename OutImageT::x_iterator cnvPtr = convolvedImage.x_at(cnvStartX, i + cnvStartY),
201  cnvEnd = cnvPtr + cnvWidth;
202  cnvPtr != cnvEnd; ++cnvPtr, ++inPtr) {
203  *cnvPtr = *inPtr;
204  }
205  }
206 }
207 
208 template <typename OutImageT, typename InImageT>
209 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage,
210  math::LinearCombinationKernel const& kernel,
211  math::ConvolutionControl const& convolutionControl) {
212  if (!kernel.isSpatiallyVarying()) {
213  // use the standard algorithm for the spatially invariant case
214  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
215  "basicConvolve for LinearCombinationKernel: spatially invariant; using brute force");
216  return convolveWithBruteForce(convolvedImage, inImage, kernel, convolutionControl.getDoNormalize());
217  } else {
218  // refactor the kernel if this is reasonable and possible;
219  // then use the standard algorithm for the spatially varying case
220  std::shared_ptr<Kernel> refKernelPtr; // possibly refactored version of kernel
221  if (static_cast<int>(kernel.getNKernelParameters()) > kernel.getNSpatialParameters()) {
222  // refactoring will speed convolution, so try it
223  refKernelPtr = kernel.refactor();
224  if (!refKernelPtr) {
225  refKernelPtr = kernel.clone();
226  }
227  } else {
228  // too few basis kernels for refactoring to be worthwhile
229  refKernelPtr = kernel.clone();
230  }
231  if (convolutionControl.getMaxInterpolationDistance() > 1) {
232  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
233  "basicConvolve for LinearCombinationKernel: using interpolation");
234  return convolveWithInterpolation(convolvedImage, inImage, *refKernelPtr, convolutionControl);
235  } else {
236  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
237  "basicConvolve for LinearCombinationKernel: maxInterpolationError < 0; using brute "
238  "force");
239  return convolveWithBruteForce(convolvedImage, inImage, *refKernelPtr,
240  convolutionControl.getDoNormalize());
241  }
242  }
243 }
244 
245 template <typename OutImageT, typename InImageT>
246 void basicConvolve(OutImageT& convolvedImage, InImageT const& inImage, math::SeparableKernel const& kernel,
247  math::ConvolutionControl const& convolutionControl) {
248  using KernelPixel = typename math::Kernel::Pixel;
249  using KernelVector = typename std::vector<KernelPixel>;
250  using KernelIterator = KernelVector::const_iterator;
251  using InXIterator = typename InImageT::const_x_iterator;
252  using InXYLocator = typename InImageT::const_xy_locator;
253  using OutXIterator = typename OutImageT::x_iterator;
254  using OutYIterator = typename OutImageT::y_iterator;
255  using OutPixel = typename OutImageT::SinglePixel;
256 
257  assertDimensionsOK(convolvedImage, inImage, kernel);
258 
259  lsst::geom::Box2I const fullBBox = inImage.getBBox(image::LOCAL);
260  lsst::geom::Box2I const goodBBox = kernel.shrinkBBox(fullBBox);
261 
262  KernelVector kernelXVec(kernel.getWidth());
263  KernelVector kernelYVec(kernel.getHeight());
264 
265  if (kernel.isSpatiallyVarying()) {
266  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
267  "SeparableKernel basicConvolve: kernel is spatially varying");
268 
269  for (int cnvY = goodBBox.getMinY(); cnvY <= goodBBox.getMaxY(); ++cnvY) {
270  double const rowPos = inImage.indexToPosition(cnvY, image::Y);
271 
272  InXYLocator inImLoc = inImage.xy_at(0, cnvY - goodBBox.getMinY());
273  OutXIterator cnvXIter = convolvedImage.row_begin(cnvY) + goodBBox.getMinX();
274  for (int cnvX = goodBBox.getMinX(); cnvX <= goodBBox.getMaxX();
275  ++cnvX, ++inImLoc.x(), ++cnvXIter) {
276  double const colPos = inImage.indexToPosition(cnvX, image::X);
277 
278  KernelPixel kSum = kernel.computeVectors(kernelXVec, kernelYVec,
279  convolutionControl.getDoNormalize(), colPos, rowPos);
280 
281  // why does this trigger warnings? It did not in the past.
282  *cnvXIter = math::convolveAtAPoint<OutImageT, InImageT>(inImLoc, kernelXVec, kernelYVec);
283  if (convolutionControl.getDoNormalize()) {
284  *cnvXIter = *cnvXIter / kSum;
285  }
286  }
287  }
288  } else {
289  // kernel is spatially invariant
290  // The basic sequence:
291  // - For each output row:
292  // - Compute x-convolved data: a kernel height's strip of input image convolved with kernel x vector
293  // - Compute one row of output by dotting each column of x-convolved data with the kernel y vector
294  // The x-convolved data is stored in a kernel-height by good-width buffer.
295  // This is circular buffer along y (to avoid shifting pixels before setting each new row);
296  // so for each new row the kernel y vector is rotated to match the order of the x-convolved data.
297 
298  LOGL_DEBUG("TRACE2.lsst.afw.math.convolve.basicConvolve",
299  "SeparableKernel basicConvolve: kernel is spatially invariant");
300 
301  kernel.computeVectors(kernelXVec, kernelYVec, convolutionControl.getDoNormalize());
302  KernelIterator const kernelXVecBegin = kernelXVec.begin();
303  KernelIterator const kernelYVecBegin = kernelYVec.begin();
304 
305  // buffer for x-convolved data
306  OutImageT buffer(lsst::geom::Extent2I(goodBBox.getWidth(), kernel.getHeight()));
307 
308  // pre-fill x-convolved data buffer with all but one row of data
309  int yInd = 0; // during initial fill bufY = inImageY
310  int const yPrefillEnd = buffer.getHeight() - 1;
311  for (; yInd < yPrefillEnd; ++yInd) {
312  OutXIterator bufXIter = buffer.x_at(0, yInd);
313  OutXIterator const bufXEnd = buffer.x_at(goodBBox.getWidth(), yInd);
314  InXIterator inXIter = inImage.x_at(0, yInd);
315  for (; bufXIter != bufXEnd; ++bufXIter, ++inXIter) {
316  *bufXIter = kernelDotProduct<OutPixel, InXIterator, KernelIterator, KernelPixel>(
317  inXIter, kernelXVecBegin, kernel.getWidth());
318  }
319  }
320 
321  // compute output pixels using the sequence described above
322  int inY = yPrefillEnd;
323  int bufY = yPrefillEnd;
324  int cnvY = goodBBox.getMinY();
325  while (true) {
326  // fill next buffer row and compute output row
327  InXIterator inXIter = inImage.x_at(0, inY);
328  OutXIterator bufXIter = buffer.x_at(0, bufY);
329  OutXIterator cnvXIter = convolvedImage.x_at(goodBBox.getMinX(), cnvY);
330  for (int bufX = 0; bufX < goodBBox.getWidth(); ++bufX, ++cnvXIter, ++bufXIter, ++inXIter) {
331  // note: bufXIter points to the row of the buffer that is being updated,
332  // whereas bufYIter points to row 0 of the buffer
333  *bufXIter = kernelDotProduct<OutPixel, InXIterator, KernelIterator, KernelPixel>(
334  inXIter, kernelXVecBegin, kernel.getWidth());
335 
336  OutYIterator bufYIter = buffer.y_at(bufX, 0);
337  *cnvXIter = kernelDotProduct<OutPixel, OutYIterator, KernelIterator, KernelPixel>(
338  bufYIter, kernelYVecBegin, kernel.getHeight());
339  }
340 
341  // test for done now, instead of the start of the loop,
342  // to avoid an unnecessary extra rotation of the kernel Y vector
343  if (cnvY >= goodBBox.getMaxY()) break;
344 
345  // update y indices, including bufY, and rotate the kernel y vector to match
346  ++inY;
347  bufY = (bufY + 1) % kernel.getHeight();
348  ++cnvY;
349  std::rotate(kernelYVec.begin(), kernelYVec.end() - 1, kernelYVec.end());
350  }
351  }
352 }
353 
354 template <typename OutImageT, typename InImageT>
355 void convolveWithBruteForce(OutImageT& convolvedImage, InImageT const& inImage, math::Kernel const& kernel,
356  math::ConvolutionControl const& convolutionControl) {
357  bool doNormalize = convolutionControl.getDoNormalize();
358 
359  using KernelPixel = typename math::Kernel::Pixel;
360  using KernelImage = image::Image<KernelPixel>;
361 
362  using KernelXIterator = typename KernelImage::const_x_iterator;
363  using KernelXYLocator = typename KernelImage::const_xy_locator;
364  using InXIterator = typename InImageT::const_x_iterator;
365  using InXYLocator = typename InImageT::const_xy_locator;
366  using OutXIterator = typename OutImageT::x_iterator;
367  using OutPixel = typename OutImageT::SinglePixel;
368 
369  assertDimensionsOK(convolvedImage, inImage, kernel);
370 
371  int const inImageWidth = inImage.getWidth();
372  int const inImageHeight = inImage.getHeight();
373  int const kWidth = kernel.getWidth();
374  int const kHeight = kernel.getHeight();
375  int const cnvWidth = inImageWidth + 1 - kernel.getWidth();
376  int const cnvHeight = inImageHeight + 1 - kernel.getHeight();
377  int const cnvStartX = kernel.getCtr().getX();
378  int const cnvStartY = kernel.getCtr().getY();
379  int const cnvEndX = cnvStartX + cnvWidth; // end index + 1
380  int const cnvEndY = cnvStartY + cnvHeight; // end index + 1
381 
382  KernelImage kernelImage(kernel.getDimensions());
383  KernelXYLocator const kernelLoc = kernelImage.xy_at(0, 0);
384 
385  if (kernel.isSpatiallyVarying()) {
386  LOGL_DEBUG("TRACE4.lsst.afw.math.convolve.convolveWithBruteForce",
387  "convolveWithBruteForce: kernel is spatially varying");
388 
389  for (int cnvY = cnvStartY; cnvY != cnvEndY; ++cnvY) {
390  double const rowPos = inImage.indexToPosition(cnvY, image::Y);
391 
392  InXYLocator inImLoc = inImage.xy_at(0, cnvY - cnvStartY);
393  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
394  for (int cnvX = cnvStartX; cnvX != cnvEndX; ++cnvX, ++inImLoc.x(), ++cnvXIter) {
395  double const colPos = inImage.indexToPosition(cnvX, image::X);
396 
397  KernelPixel kSum = kernel.computeImage(kernelImage, false, colPos, rowPos);
398  *cnvXIter = math::convolveAtAPoint<OutImageT, InImageT>(inImLoc, kernelLoc, kWidth, kHeight);
399  if (doNormalize) {
400  *cnvXIter = *cnvXIter / kSum;
401  }
402  }
403  }
404  } else {
405  LOGL_DEBUG("TRACE4.lsst.afw.math.convolve.convolveWithBruteForce",
406  "convolveWithBruteForce: kernel is spatially invariant");
407 
408  (void)kernel.computeImage(kernelImage, doNormalize);
409 
410  for (int inStartY = 0, cnvY = cnvStartY; inStartY < cnvHeight; ++inStartY, ++cnvY) {
411  KernelXIterator kernelXIter = kernelImage.x_at(0, 0);
412  InXIterator inXIter = inImage.x_at(0, inStartY);
413  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
414  for (int x = 0; x < cnvWidth; ++x, ++cnvXIter, ++inXIter) {
415  *cnvXIter = kernelDotProduct<OutPixel, InXIterator, KernelXIterator, KernelPixel>(
416  inXIter, kernelXIter, kWidth);
417  }
418  for (int kernelY = 1, inY = inStartY + 1; kernelY < kHeight; ++inY, ++kernelY) {
419  KernelXIterator kernelXIter = kernelImage.x_at(0, kernelY);
420  InXIterator inXIter = inImage.x_at(0, inY);
421  OutXIterator cnvXIter = convolvedImage.x_at(cnvStartX, cnvY);
422  for (int x = 0; x < cnvWidth; ++x, ++cnvXIter, ++inXIter) {
423  *cnvXIter += kernelDotProduct<OutPixel, InXIterator, KernelXIterator, KernelPixel>(
424  inXIter, kernelXIter, kWidth);
425  }
426  }
427  }
428  }
429 }
430 
431 /*
432  * Explicit instantiation
433  */
435 #define IMAGE(PIXTYPE) image::Image<PIXTYPE>
436 #define MASKEDIMAGE(PIXTYPE) image::MaskedImage<PIXTYPE, image::MaskPixel, image::VariancePixel>
437 #define NL /* */
438 // Instantiate Image or MaskedImage versions
439 #define INSTANTIATE_IM_OR_MI(IMGMACRO, OUTPIXTYPE, INPIXTYPE) \
440  template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, math::Kernel const&, \
441  math::ConvolutionControl const&); \
442  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
443  math::DeltaFunctionKernel const&, math::ConvolutionControl const&); \
444  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
445  math::LinearCombinationKernel const&, math::ConvolutionControl const&); \
446  NL template void basicConvolve(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
447  math::SeparableKernel const&, math::ConvolutionControl const&); \
448  NL template void convolveWithBruteForce(IMGMACRO(OUTPIXTYPE)&, IMGMACRO(INPIXTYPE) const &, \
449  math::Kernel const&, math::ConvolutionControl const&);
450 // Instantiate both Image and MaskedImage versions
451 #define INSTANTIATE(OUTPIXTYPE, INPIXTYPE) \
452  INSTANTIATE_IM_OR_MI(IMAGE, OUTPIXTYPE, INPIXTYPE) \
453  INSTANTIATE_IM_OR_MI(MASKEDIMAGE, OUTPIXTYPE, INPIXTYPE)
454 
455 INSTANTIATE(double, double)
456 INSTANTIATE(double, float)
457 INSTANTIATE(double, int)
458 INSTANTIATE(double, std::uint16_t)
459 INSTANTIATE(float, float)
460 INSTANTIATE(float, int)
462 INSTANTIATE(int, int)
465 } // namespace detail
466 } // namespace math
467 } // namespace afw
468 } // namespace lsst
double x
#define IS_INSTANCE(A, B)
Definition: Convolve.h:40
#define INSTANTIATE(FROMSYS, TOSYS)
Definition: Detector.cc:484
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
LSST DM logging module built on log4cxx.
#define LOGL_DEBUG(logger, message...)
Log a debug-level message using a varargs/printf style interface.
Definition: Log.h:515
std::ostream * os
Definition: Schema.cc:557
A class to represent a 2-dimensional array of pixels.
Definition: Image.h:51
Parameters to control convolution.
Definition: ConvolveImage.h:50
A kernel that has only one non-zero pixel (of value 1)
Definition: Kernel.h:643
lsst::geom::Point2I getPixel() const
Definition: Kernel.h:670
Kernels are used for convolution with MaskedImages and (eventually) Images.
Definition: Kernel.h:110
lsst::geom::Extent2I const getDimensions() const
Return the Kernel's dimensions (width, height)
Definition: Kernel.h:212
int getHeight() const
Return the Kernel's height.
Definition: Kernel.h:229
lsst::geom::Point2I getCtr() const
Return index of kernel's center.
Definition: Kernel.h:234
unsigned int getNKernelParameters() const
Return the number of kernel parameters (0 if none)
Definition: Kernel.h:246
int getNSpatialParameters() const
Return the number of spatial parameters (0 if not spatially varying)
Definition: Kernel.h:251
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:173
int getWidth() const
Return the Kernel's width.
Definition: Kernel.h:224
bool isSpatiallyVarying() const
Return true iff the kernel is spatially varying (has a spatial function)
Definition: Kernel.h:333
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:76
A kernel that is a linear combination of fixed basis kernels.
Definition: Kernel.h:704
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...
std::shared_ptr< Kernel > clone() const override
Return a pointer to a deep copy of this kernel.
A kernel described by a pair of functions: func(x, y) = colFunc(x) * rowFunc(y)
Definition: Kernel.h:860
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)
An integer coordinate rectangle.
Definition: Box.h:55
int getMinY() const noexcept
Definition: Box.h:158
int getMinX() const noexcept
Definition: Box.h:157
int getWidth() const noexcept
Definition: Box.h:187
int getMaxX() const noexcept
Definition: Box.h:161
int getMaxY() const noexcept
Definition: Box.h:162
Reports invalid arguments.
Definition: Runtime.h:66
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.
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.
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.
A base class for image defects.
T rotate(T... args)