LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
BasisLists.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
11 #include <cmath>
12 #include <limits>
13 
14 #include "boost/timer.hpp"
15 
17 #include "lsst/pex/policy/Policy.h"
18 #include "lsst/pex/logging/Trace.h"
19 #include "lsst/afw/geom.h"
20 #include "lsst/afw/image.h"
21 #include "lsst/afw/math.h"
23 
24 namespace pexExcept = lsst::pex::exceptions;
25 namespace afwGeom = lsst::afw::geom;
26 namespace afwImage = lsst::afw::image;
27 namespace afwMath = lsst::afw::math;
28 namespace pexLogging = lsst::pex::logging;
29 
30 namespace lsst {
31 namespace ip {
32 namespace diffim {
33 
49  int width,
50  int height
51  ) {
52  if ((width < 1) || (height < 1)) {
53  throw LSST_EXCEPT(pexExcept::Exception, "nRows and nCols must be positive");
54  }
55  const int signedWidth = static_cast<int>(width);
56  const int signedHeight = static_cast<int>(height);
57  afwMath::KernelList kernelBasisList;
58  for (int row = 0; row < signedHeight; ++row) {
59  for (int col = 0; col < signedWidth; ++col) {
60  boost::shared_ptr<afwMath::Kernel>
61  kernelPtr(new afwMath::DeltaFunctionKernel(width, height, afwGeom::Point2I(col,row)));
62  kernelBasisList.push_back(kernelPtr);
63  }
64  }
65  return kernelBasisList;
66  }
67 
80  int halfWidth,
81  int nGauss,
82  std::vector<double> const &sigGauss,
83  std::vector<int> const &degGauss
84  ) {
85  typedef afwMath::Kernel::Pixel Pixel;
86  typedef afwImage::Image<Pixel> Image;
87 
88  if (halfWidth < 1) {
89  throw LSST_EXCEPT(pexExcept::Exception, "halfWidth must be positive");
90  }
91  if (nGauss != static_cast<int>(sigGauss.size())) {
92  throw LSST_EXCEPT(pexExcept::Exception, "sigGauss does not have enough entries");
93  }
94  if (nGauss != static_cast<int>(degGauss.size())) {
95  throw LSST_EXCEPT(pexExcept::Exception, "degGauss does not have enough entries");
96  }
97  int fullWidth = 2 * halfWidth + 1;
98  Image image(afwGeom::Extent2I(fullWidth, fullWidth));
99 
100  afwMath::KernelList kernelBasisList;
101  for (int i = 0; i < nGauss; i++) {
102  /*
103  sigma = FWHM / ( 2 * sqrt(2 * ln(2)) )
104  */
105  double sig = sigGauss[i];
106  int deg = degGauss[i];
107 
108  pexLogging::TTrace<2>("lsst.ip.diffim.BasisLists.makeAlardLuptonBasisList",
109  "Gaussian %d : sigma %.2f degree %d", i, sig, deg);
110 
111  afwMath::GaussianFunction2<Pixel> gaussian(sig, sig);
112  afwMath::AnalyticKernel kernel(fullWidth, fullWidth, gaussian);
113  afwMath::PolynomialFunction2<Pixel> polynomial(deg);
114 
115  for (int j = 0, n = 0; j <= deg; j++) {
116  for (int k = 0; k <= (deg - j); k++, n++) {
117  /* for 0th order term, skip polynomial */
118  (void)kernel.computeImage(image, true);
119  if (n == 0) {
120  boost::shared_ptr<afwMath::Kernel>
121  kernelPtr(new afwMath::FixedKernel(image));
122  kernelBasisList.push_back(kernelPtr);
123  continue;
124  }
125 
126  /* gaussian to be modified by this term in the polynomial */
127  polynomial.setParameter(n, 1.);
128  (void)kernel.computeImage(image, true);
129  for (int y = 0, v = -halfWidth; y < image.getHeight(); y++, v++) {
130  int u = -halfWidth;
131  for (Image::xy_locator ptr = image.xy_at(0, y),
132  end = image.xy_at(image.getWidth(), y);
133  ptr != end; ++ptr.x(), u++) {
134  /* Evaluate from -1 to 1 */
135  *ptr = *ptr * polynomial(u/static_cast<double>(halfWidth),
136  v/static_cast<double>(halfWidth));
137  }
138  }
139  boost::shared_ptr<afwMath::Kernel>
140  kernelPtr(new afwMath::FixedKernel(image));
141  kernelBasisList.push_back(kernelPtr);
142  polynomial.setParameter(n, 0.);
143  }
144  }
145  }
146  return renormalizeKernelList(kernelBasisList);
147  }
148 
149 
150  boost::shared_ptr<Eigen::MatrixXd>
153  ) {
154 
155  /* NOTES
156  *
157  * The 3-point first derivative central difference (Laplacian) yields
158  * diagonal stripes and should not be used.
159 
160  coeffs[0][1] = -1. / 2.;
161  coeffs[1][0] = -1. / 2.;
162  coeffs[1][2] = 1. / 2.;
163  coeffs[2][1] = 1. / 2.;
164 
165  * The 5-point second derivative central difference (Laplacian) looks great.
166  * For smaller lambdas you need a higher border penalty. In general, if you
167  * decrease the regularization strength you should increase the border
168  * penalty to avoid noise in the border pixels. This has been used to value
169  * 10 and things still look OK.
170 
171  * The 9-point second derivative central difference (Laplacian with off
172  * diagonal terms) also looks great. Seems to have slightly higher-valued
173  * border pixels so make boundary larger if using this. E.g. 1.5.
174 
175  * The forward finite difference, first derivative term works good
176  *
177  * The forward finite difference, second derivative term is slightly banded
178  * in LLC
179  *
180  * The forward finite difference, third derivative term is highly banded in
181  * LLC and should not be used.
182  *
183  */
184 
185  std::string regularizationType = policy.getString("regularizationType");
186  int width = policy.getInt("kernelSize");
187  int height = policy.getInt("kernelSize");
188  float borderPenalty = policy.getDouble("regularizationBorderPenalty");
189  bool fitForBackground = policy.getBool("fitForBackground");
190 
191  boost::shared_ptr<Eigen::MatrixXd> bMat;
192  if (regularizationType == "centralDifference") {
193  int stencil = policy.getInt("centralRegularizationStencil");
194  bMat = makeCentralDifferenceMatrix(width, height, stencil, borderPenalty, fitForBackground);
195  }
196  else if (regularizationType == "forwardDifference") {
197  std::vector<int> orders = policy.getIntArray("forwardRegularizationOrders");
198  bMat = makeForwardDifferenceMatrix(width, height, orders, borderPenalty, fitForBackground);
199  }
200  else {
201  throw LSST_EXCEPT(pexExcept::Exception, "regularizationType not recognized");
202  }
203 
204  boost::shared_ptr<Eigen::MatrixXd> hMat (new Eigen::MatrixXd((*bMat).transpose() * (*bMat)));
205  return hMat;
206  }
207 
211  boost::shared_ptr<Eigen::MatrixXd>
213  int width,
214  int height,
215  int stencil,
216  float borderPenalty,
217  bool fitForBackground
218  ) {
219 
220  /* 5- or 9-point stencil to approximate the Laplacian; i.e. this is a second
221  * order central finite difference.
222  *
223  * 5-point approximation ignores off-diagonal elements, and is
224  *
225  * f(x-1,y) + f(x+1,y) + f(x,y-1) + f(x,y+1) - 4 f(x,y)
226  *
227  * 0 1 0
228  * 1 -4 1
229  * 0 1 0
230  *
231  *
232  * 9-point approximation includes diagonals and is
233  *
234  * 1 4 1
235  * 4 -20 4 / 6
236  * 1 4 1
237  *
238  */
239 
240  if (borderPenalty < 0)
241  throw LSST_EXCEPT(pexExcept::Exception, "Only border penalty of >= 0 allowed");
242 
243  std::vector<std::vector<float> >
244  coeffs(3, std::vector<float>(3,0));
245 
246  if (stencil == 5) {
247  /* http://www.physics.arizona.edu/~restrepo/475B/Notes/source/node51.html
248  *
249  * This is equivalent to a second order central difference summed along
250  * the x-axis, and then summed along the y-axis.
251  *
252  * http://www.holoborodko.com/pavel/?page_id=239
253  *
254  */
255  coeffs[0][1] = 1.;
256  coeffs[1][0] = 1.;
257  coeffs[1][1] = -4.;
258  coeffs[1][2] = 1.;
259  coeffs[2][1] = 1.;
260  }
261  else if (stencil == 9) {
262  /* http://www.physics.arizona.edu/~restrepo/475B/Notes/source/node52.html */
263  coeffs[0][0] = 1. / 6.;
264  coeffs[0][1] = 4. / 6.;
265  coeffs[0][2] = 1. / 6.;
266  coeffs[1][0] = 4. / 6.;
267  coeffs[1][1] = -20. / 6.;
268  coeffs[1][2] = 4. / 6.;
269  coeffs[2][0] = 1. / 6.;
270  coeffs[2][1] = 4. / 6.;
271  coeffs[2][2] = 1. / 6.;
272  }
273  else {
274  throw LSST_EXCEPT(pexExcept::Exception, "Only 5- or 9-point Laplacian stencils allowed");
275  }
276 
277  int nBgTerms = fitForBackground ? 1 : 0;
278  Eigen::MatrixXd bMat = Eigen::MatrixXd::Zero(width * height + nBgTerms, width * height + nBgTerms);
279 
280  for (int i = 0; i < width*height; i++) {
281  int const x0 = i % width; // the x coord in the kernel image
282  int const y0 = i / width; // the y coord in the kernel image
283  int const distX = width - x0 - 1; // distance from edge of image
284  int const distY = height - y0 - 1; // distance from edge of image
285 
286  if ( (x0 > 0) && (y0 > 0) && (distX > 0) && (distY > 0) ) {
287  for (int dx = -1; dx < 2; dx += 1) {
288  for (int dy = -1; dy < 2; dy += 1) {
289  bMat(i, i + dx + dy * width) += coeffs[dx+1][dy+1];
290  }
291  }
292  }
293  else {
294  bMat(i, i) = borderPenalty;
295  }
296  }
297 
298  if (fitForBackground) {
299  /* Last row / col should have no regularization since its the background term */
300  if (bMat.col(width*height).sum() != 0.) {
301  throw LSST_EXCEPT(pexExcept::Exception, "Error 1 in regularization matrix");
302  }
303  if (bMat.row(width*height).sum() != 0.) {
304  throw LSST_EXCEPT(pexExcept::Exception, "Error 2 in regularization matrix");
305  }
306  }
307 
308  boost::shared_ptr<Eigen::MatrixXd> bMatPtr (new Eigen::MatrixXd(bMat));
309  return bMatPtr;
310  }
311 
315  boost::shared_ptr<Eigen::MatrixXd>
317  int width,
318  int height,
319  std::vector<int> const& orders,
320  float borderPenalty,
321  bool fitForBackground
322  ) {
323 
324  /*
325  Instead of Taylor expanding the forward difference approximation of
326  derivatives (N.R. section 18.5) lets just hard code in the expansion of
327  the 1st through 3rd derivatives, which will try and enforce smoothness of
328  0 through 2nd derivatives.
329 
330  A property of the basic "finite difference regularization" is that their
331  rows (column multipliers) sum to 0.
332 
333  We taper the order of the constraint as you get close to the boundary.
334  The actual perimeter values are left unconstrained but we might want to
335  consider giving these the same penalties as the other border pixels,
336  which is +1 (since the value gets squared).
337 
338  */
339 
340  if (borderPenalty < 0)
341  throw LSST_EXCEPT(pexExcept::Exception, "Only border penalty of >= 0 allowed");
342 
343  std::vector<std::vector<float> >
344  coeffs(4, std::vector<float>(4,0));
345 
346  // penalize border? this is along the diagonal and gets squared, so sign does not matter
347  coeffs[0][0] = borderPenalty;
348 
349  // penalize zeroth derivative
350  coeffs[1][0] = -1.;
351  coeffs[1][1] = +1.;
352 
353  // penalize first derivative
354  coeffs[2][0] = -1.;
355  coeffs[2][1] = +2.;
356  coeffs[2][2] = -1.;
357 
358  // penalize second derivative
359  coeffs[3][0] = -1.;
360  coeffs[3][1] = +3.;
361  coeffs[3][2] = -3.;
362  coeffs[3][3] = +1.;
363 
364  int nBgTerms = fitForBackground ? 1 : 0;
365  Eigen::MatrixXd bTot = Eigen::MatrixXd::Zero(width * height + nBgTerms, width * height + nBgTerms);
366 
367  std::vector<int>::const_iterator order;
368  for (order = orders.begin(); order != orders.end(); order++) {
369  if ((*order < 1) || (*order > 3))
370  throw LSST_EXCEPT(pexExcept::Exception, "Only orders 1..3 allowed");
371 
372  Eigen::MatrixXd bMatX = Eigen::MatrixXd::Zero(width * height + nBgTerms,
373  width * height + nBgTerms);
374  Eigen::MatrixXd bMatY = Eigen::MatrixXd::Zero(width * height + nBgTerms,
375  width * height + nBgTerms);
376 
377  for (int i = 0; i < width*height; i++) {
378  int const x0 = i % width; // the x coord in the kernel image
379  int const y0 = i / width; // the y coord in the kernel image
380 
381  int distX = width - x0 - 1; // distance from edge of image
382  int orderToUseX = std::min(distX, *order);
383  for (int j = 0; j < orderToUseX+1; j++) {
384  bMatX(i, i + j) = coeffs[orderToUseX][j];
385  }
386 
387  int distY = height - y0 - 1; // distance from edge of image
388  int orderToUseY = std::min(distY, *order);
389  for (int j = 0; j < orderToUseY+1; j++) {
390  bMatY(i, i + j * width) = coeffs[orderToUseY][j];
391  }
392  }
393  bTot += bMatX;
394  bTot += bMatY;
395  }
396 
397  if (fitForBackground) {
398  /* Last row / col should have no regularization since its the background term */
399  if (bTot.col(width*height).sum() != 0.) {
400  throw LSST_EXCEPT(pexExcept::Exception, "Error in regularization matrix");
401  }
402  if (bTot.row(width*height).sum() != 0.) {
403  throw LSST_EXCEPT(pexExcept::Exception, "Error in regularization matrix");
404  }
405  }
406 
407  boost::shared_ptr<Eigen::MatrixXd> bMatPtr (new Eigen::MatrixXd(bTot));
408  return bMatPtr;
409  }
410 
411 
421  lsst::afw::math::KernelList const &kernelListIn
422  ) {
423  typedef afwMath::Kernel::Pixel Pixel;
424  typedef afwImage::Image<Pixel> Image;
425  double kSum;
426 
427  /*
428 
429  We want all the bases except for the first to sum to 0.0. This allows
430  us to achieve kernel flux conservation (Ksum) across the image since all
431  the power will be in the first term, which will not vary spatially.
432 
433  K(x,y) = Ksum * B_0 + Sum_i : a(x,y) * B_i
434 
435  To do this, normalize all Kernels to sum = 1. and subtract B_0 from all
436  subsequent kenrels.
437 
438  To get an idea of the relative contribution of each of these basis
439  functions later on down the line, lets also normalize them such that
440 
441  Sum(B_i) == 0.0 *and*
442  B_i * B_i == 1.0
443 
444  For completeness
445 
446  Sum(B_0) == 1.0
447  B_0 * B_0 != 1.0
448 
449  */
450  afwMath::KernelList kernelListOut;
451  if (kernelListIn.size() == 0) {
452  return kernelListOut;
453  }
454 
455  Image image0(kernelListIn[0]->getDimensions());
456  for (unsigned int i = 0; i < kernelListIn.size(); i++) {
457  if (i == 0) {
458  /* Make sure that it is normalized to kSum 1. */
459  (void)kernelListIn[i]->computeImage(image0, true);
460  boost::shared_ptr<afwMath::Kernel>
461  kernelPtr(new afwMath::FixedKernel(image0));
462  kernelListOut.push_back(kernelPtr);
463 
464  continue;
465  }
466 
467  /* Don't normalize here */
468  Image image(kernelListIn[i]->getDimensions());
469  (void)kernelListIn[i]->computeImage(image, false);
470  /* image.writeFits(str(boost::format("in_k%d.fits") % i)); */
471 
472  /* Check the kernel sum; if its close to zero don't do anything */
473  kSum = 0.;
474  for (int y = 0; y < image.getHeight(); y++) {
475  for (Image::xy_locator ptr = image.xy_at(0, y), end = image.xy_at(image.getWidth(), y);
476  ptr != end; ++ptr.x()) {
477  kSum += *ptr;
478  }
479  }
480 
481  /* std::numeric_limits<float>::epsilon() ~ e-7
482  std::numeric_limits<double>::epsilon() ~ e-16
483 
484  If we end up with 2e-16 kernel sum, this still blows up the kernel values.
485  Even though the kernels are double, use the float limits instead
486  */
487  if (fabs(kSum) > std::numeric_limits<float>::epsilon()) {
488  image /= kSum;
489  image -= image0;
490  }
491 
492  /* Finally, rescale such that the inner product is 1 */
493  kSum = 0.;
494  for (int y = 0; y < image.getHeight(); y++) {
495  for (Image::xy_locator ptr = image.xy_at(0, y), end = image.xy_at(image.getWidth(), y);
496  ptr != end; ++ptr.x()) {
497  kSum += *ptr * *ptr;
498  }
499  }
500  image /= std::sqrt(kSum);
501  /* image.writeFits(str(boost::format("out_k%d.fits") % i)); */
502 
503  boost::shared_ptr<afwMath::Kernel>
504  kernelPtr(new afwMath::FixedKernel(image));
505  kernelListOut.push_back(kernelPtr);
506  }
507  return kernelListOut;
508  }
509 
510 
511 
515  boost::shared_ptr<Eigen::MatrixXd>
517  unsigned int width,
518  unsigned int height,
519  unsigned int order,
520  unsigned int boundary_style, // 0 = unwrapped, 1 = wrapped, 2 = order-tapered
521  unsigned int difference_style, // 0 = forward, 1 = central
522  bool printB // a debug flag ... remove when done.
523  ) {
524 
525  if (order > 2) throw LSST_EXCEPT(pexExcept::Exception, "Only orders 0..2 allowed");
526 
527  if (boundary_style > 2) {
528  throw LSST_EXCEPT(pexExcept::Exception, "Boundary styles 0..2 defined");
529  }
530  if (difference_style > 1) {
531  throw LSST_EXCEPT(pexExcept::Exception, "Only forward(0), and central(1) difference types defined");
532  }
533 
534  /* what works, and what doesn't */
535  // == good job ==
536  // - order 0, wrapped, forward
537  // - order 1, wrapped or tapered, central or forward
538  // - order 2, wrapped or tapered, forward
539  // == bad job (usually diagongal stripes) ==
540  // - all others
541 
542 
543  /*
544  Instead of Taylor expanding the forward difference approximation of
545  derivatives (N.R. section 18.5) lets just hard code in the expansion of
546  the 1st through 3rd derivatives, which will try and enforce smoothness of
547  0 through 2nd derivatives.
548 
549  A property of the basic "finite difference regularization" is that their
550  rows (column multipliers) sum to 0.
551 
552  Another consideration is to use *multiple* finite difference operators as
553  a constraint.
554 
555  */
556 
557 
558  // ===========================================================================
559  // Get the coeffs for the finite differencing
560  // note: The coeffs are stored 2D although they are essentially 1D entities.
561  // The 2d design was chosen to allow cross-terms to be included,
562  // though they are not yet implemented.
563  //
564  std::vector<std::vector<std::vector<float> > >
565  coeffs(3, std::vector<std::vector<float> >(5, std::vector<float>(5,0)));
566  unsigned int x_cen = 0, y_cen = 0; // center of reqested order coeffs
567  unsigned int x_cen1 = 0, y_cen1 = 0; // center of order 1 coeffs
568  unsigned int x_cen2 = 0, y_cen2 = 0; // center of order 2 coeffs
569  unsigned int x_size = 0, y_size = 0;
570 
571  // forward difference coefficients
572  if (difference_style == 0) {
573 
574  y_cen = x_cen = 0;
575  x_cen1 = y_cen1 = 0;
576  x_cen2 = y_cen2 = 0;
577 
578  x_size = y_size = order + 2;
579 
580  // default forward difference suggested in NR chap 18
581  // 0th order
582  coeffs[0][0][0] = -2;
583  coeffs[0][0][1] = 1;
584  coeffs[0][1][0] = 1;
585  coeffs[0][1][1] = 0;
586 
587  // 1st 2
588  coeffs[1][0][0] = -2;
589  coeffs[1][0][1] = 2;
590  coeffs[1][0][2] = -1;
591  coeffs[1][1][0] = 2;
592  coeffs[1][1][1] = 0;
593  coeffs[1][1][2] = 0;
594  coeffs[1][2][0] = -1;
595  coeffs[1][2][1] = 0;
596  coeffs[1][2][2] = 0;
597 
598  // 2nd 2
599  coeffs[2][0][0] = -2;
600  coeffs[2][0][1] = 3;
601  coeffs[2][0][2] = -3;
602  coeffs[2][0][3] = 1;
603  coeffs[2][1][0] = 3;
604  coeffs[2][1][1] = 0;
605  coeffs[2][1][2] = 0;
606  coeffs[2][1][3] = 0;
607  coeffs[2][2][0] = -3;
608  coeffs[2][2][1] = 0;
609  coeffs[2][2][2] = 0;
610  coeffs[2][2][3] = 0;
611  coeffs[2][3][0] = 1;
612  coeffs[2][3][1] = 0;
613  coeffs[2][3][2] = 0;
614  coeffs[2][3][3] = 0;
615  }
616 
617  // central difference coefficients
618  if (difference_style == 1) {
619 
620  // this is asymmetric and produces diagonal banding in the kernel
621  // from: http://www.holoborodko.com/pavel/?page_id=239
622  if (order == 0) {
623  y_cen = x_cen = 1;
624  x_size = y_size = 3;
625  }
626  coeffs[0][0][0] = 0;
627  coeffs[0][0][1] = -1;
628  coeffs[0][0][2] = 0;
629 
630  coeffs[0][1][0] = -1;
631  coeffs[0][1][1] = 0;
632  coeffs[0][1][2] = 1;
633 
634  coeffs[0][2][0] = 0;
635  coeffs[0][2][1] = 1;
636  coeffs[0][2][2] = 0;
637 
638 
639  // this works well and is largely the same as order=1 forward-diff.
640  // from: http://www.holoborodko.com/pavel/?page_id=239
641  if (order == 1) {
642  y_cen = x_cen = 1;
643  x_size = y_size = 3;
644  }
645  y_cen1 = x_cen1 = 1;
646  coeffs[1][0][0] = 0;
647  coeffs[1][0][1] = 1;
648  coeffs[1][0][2] = 0;
649 
650  coeffs[1][1][0] = 1;
651  coeffs[1][1][1] = -4;
652  coeffs[1][1][2] = 1;
653 
654  coeffs[1][2][0] = 0;
655  coeffs[1][2][1] = 1;
656  coeffs[1][2][2] = 0;
657 
658 
659  // asymmetric and produces diagonal banding in the kernel
660  // from http://www.holoborodko.com/pavel/?page_id=239
661  if (order == 2) {
662  y_cen = x_cen = 2;
663  x_size = y_size = 5;
664  }
665  y_cen2 = x_cen2 = 2;
666  coeffs[2][0][0] = 0;
667  coeffs[2][0][1] = 0;
668  coeffs[2][0][2] = -1;
669  coeffs[2][0][3] = 0;
670  coeffs[2][0][4] = 0;
671 
672  coeffs[2][1][0] = 0;
673  coeffs[2][1][1] = 0;
674  coeffs[2][1][2] = 2;
675  coeffs[2][1][3] = 0;
676  coeffs[2][1][4] = 0;
677 
678  coeffs[2][2][0] = -1;
679  coeffs[2][2][1] = 2;
680  coeffs[2][2][2] = 0;
681  coeffs[2][2][3] = -2;
682  coeffs[2][2][4] = 1;
683 
684  coeffs[2][3][0] = 0;
685  coeffs[2][3][1] = 0;
686  coeffs[2][3][2] = -2;
687  coeffs[2][3][3] = 0;
688  coeffs[2][3][4] = 0;
689 
690  coeffs[2][4][0] = 0;
691  coeffs[2][4][1] = 0;
692  coeffs[2][4][2] = 1;
693  coeffs[2][4][3] = 0;
694  coeffs[2][4][4] = 0;
695  }
696 
697 
698  /* Note we have to add 1 extra (empty) term here because of the differential
699  * background fitting */
700  Eigen::MatrixXd bMat = Eigen::MatrixXd::Zero(width*height+1, width*height+1);
701 
702  /* Forward difference approximation */
703  for (unsigned int i = 0; i < width*height; i++) {
704 
705  unsigned int const x0 = i % width; // the x coord in the kernel image
706  unsigned int const y0 = i / width; // the y coord in the kernel image
707 
708  unsigned int x_edge_distance = (x0 > (width - x0 - 1)) ? width - x0 - 1 : x0;
709  unsigned int y_edge_distance = (y0 > (height - y0 - 1)) ? height - y0 - 1 : y0;
710  unsigned int edge_distance = (x_edge_distance < y_edge_distance) ? x_edge_distance :
711  y_edge_distance;
712 
713  for (unsigned int dx = 0; dx < x_size; dx++) {
714  for (unsigned int dy = 0; dy < y_size; dy++) {
715 
716  // determine where to put this coeff
717 
718  // handle the boundary condition
719  // note: adding width and height in the sum prevents negatives
720  unsigned int x = 0;
721  unsigned int y = 0;
722  double this_coeff = 0;
723 
724  // no-wrapping at edges
725  if (boundary_style == 0) {
726  x = x0 + dx - x_cen;
727  y = y0 + dy - y_cen;
728  if (y > height - 1 || x > width - 1) {
729  continue;
730  }
731  this_coeff = coeffs[order][dx][dy];
732 
733  // wrapping at edges
734  } else if (boundary_style == 1) {
735  x = (width + x0 + dx - x_cen) % width;
736  y = (height + y0 + dy - y_cen) % height;
737  this_coeff = coeffs[order][dx][dy];
738 
739  // order tapering to the edge (just clone wrapping for now)
740  // - use the lowest order possible
741  } else if (boundary_style == 2) {
742 
743  // edge rows and columns ... set to constant
744  if (edge_distance == 0) {
745  x = x0;
746  y = y0;
747  this_coeff = 1;
748  }
749  // in one from edge, use 1st order
750  else if (edge_distance == 1 && order > 0) {
751  x = (width + x0 + dx - x_cen1) % width;
752  y = (height + y0 + dy - y_cen1) % height;
753  if ((dx < 3) && (dy < 3)) { this_coeff = coeffs[1][dx][dy]; }
754  }
755  // in two from edge, use 2st order if order > 1
756  else if (edge_distance == 2 && order > 1){
757  x = (width + x0 + dx - x_cen2) % width;
758  y = (height + y0 + dy - y_cen2) % height;
759  if ((dx < 5) && (dy < 5)) { this_coeff = coeffs[2][dx][dy]; }
760  }
761  // if we're somewhere in the middle
762  else if (edge_distance > order) {
763  x = (width + x0 + dx - x_cen) % width;
764  y = (height + y0 + dy - y_cen) % height;
765  this_coeff = coeffs[order][dx][dy];
766  }
767 
768  }
769 
770  bMat(i, y*width + x) = this_coeff;
771 
772  }
773 
774  }
775 
776  }
777 
778  if (printB) {
779  std::cout << bMat << std::endl;
780  }
781 
782  boost::shared_ptr<Eigen::MatrixXd> hMat (new Eigen::MatrixXd(bMat.transpose() * bMat));
783  return hMat;
784  }
785 
786 }}} // end of namespace lsst::ip::diffim
int y
An include file to include the public header files for lsst::afw::math.
int getInt(const std::string &name) const
Definition: Policy.h:603
An include file to include the header files for lsst::afw::geom.
lsst::afw::math::KernelList renormalizeKernelList(lsst::afw::math::KernelList const &kernelListIn)
Rescale an input set of kernels.
Definition: BasisLists.cc:420
bool getBool(const std::string &name) const
Definition: Policy.h:589
a container for holding hierarchical configuration data in memory.
Definition: Policy.h:169
2-dimensional polynomial function with cross terms
const std::string getString(const std::string &name) const
Definition: Policy.h:631
definition of the Trace messaging facilities
int const x0
Definition: saturated.cc:45
boost::shared_ptr< Eigen::MatrixXd > makeCentralDifferenceMatrix(int width, int height, int stencil, float borderPenalty, bool fitForBackground)
Generate regularization matrix for delta function kernels.
Definition: BasisLists.cc:212
table::Key< table::Array< Kernel::Pixel > > image
Definition: FixedKernel.cc:117
An include file to include the header files for lsst::afw::image.
boost::shared_ptr< Eigen::MatrixXd > makeRegularizationMatrix(lsst::pex::policy::Policy policy)
Build a regularization matrix for Delta function kernels.
Definition: BasisLists.cc:151
void ImageT ImageT int float saturatedPixelValue int const width
Definition: saturated.cc:44
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.
boost::shared_ptr< Eigen::MatrixXd > makeFiniteDifferenceRegularizationDeprecated(unsigned int width, unsigned int height, unsigned int order, unsigned int boundary_style, unsigned int difference_style, bool printB)
Generate regularization matrix for delta function kernels.
Definition: BasisLists.cc:516
double x
void ImageT ImageT int float saturatedPixelValue int const height
Definition: saturated.cc:44
int row
Definition: CR.cc:153
#define LSST_EXCEPT(type,...)
Definition: Exception.h:46
double getDouble(const std::string &name) const
Definition: Policy.h:617
lsst::afw::math::KernelList makeDeltaFunctionBasisList(int width, int height)
Generate a basis set of delta function Kernels.
Definition: BasisLists.cc:48
void setParameter(unsigned int ind, double newValue)
Set one function parameter without range checking.
Definition: Function.h:161
A kernel described by a function.
Definition: Kernel.h:628
int const y0
Definition: saturated.cc:45
Subroutines associated with generating, normalising, and regularising Basis functions.
A kernel that has only one non-zero pixel (of value 1)
Definition: Kernel.h:744
int col
Definition: CR.cc:152
A kernel created from an Image.
Definition: Kernel.h:551
boost::shared_ptr< Eigen::MatrixXd > makeForwardDifferenceMatrix(int width, int height, std::vector< int > const &orders, float borderPenalty, bool fitForBackground)
Generate regularization matrix for delta function kernels.
Definition: BasisLists.cc:316
std::vector< boost::shared_ptr< Kernel > > KernelList
Definition: Kernel.h:542
IntArray getIntArray(const std::string &name) const
Definition: Policy.h:1031
lsst::afw::math::KernelList makeAlardLuptonBasisList(int halfWidth, int nGauss, std::vector< double > const &sigGauss, std::vector< int > const &degGauss)
Generate an Alard-Lupton basis set of Kernels.
Definition: BasisLists.cc:79