lsst.meas.base  14.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Macros Groups Pages
SdssCentroid.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 /*
3  * LSST Data Management System
4  * Copyright 2008-2013 LSST Corporation.
5  *
6  * This product includes software developed by the
7  * LSST Project (http://www.lsst.org/).
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the LSST License Statement and
20  * the GNU General Public License along with this program. If not,
21  * see <http://www.lsstcorp.org/LegalNotices/>.
22  */
23 #include <iostream>
24 #include <cmath>
25 #include <numeric>
26 #include "ndarray/eigen.h"
27 #include "lsst/afw/detection/Psf.h"
28 #include "lsst/pex/exceptions.h"
29 #include "lsst/afw/math/ConvolveImage.h"
30 #include "lsst/afw/math/offsetImage.h"
31 #include "lsst/afw/geom/Angle.h"
32 #include "lsst/afw/table/Source.h"
34 
35 
36 namespace lsst { namespace meas { namespace base {
37 namespace {
38 FlagDefinitionList flagDefinitions;
39 } // end anonymous
40 
41 FlagDefinition const SdssCentroidAlgorithm::FAILURE = flagDefinitions.addFailureFlag();
42 FlagDefinition const SdssCentroidAlgorithm::EDGE = flagDefinitions.add("flag_edge", "Object too close to edge");
43 FlagDefinition const SdssCentroidAlgorithm::NO_SECOND_DERIVATIVE = flagDefinitions.add("flag_noSecondDerivative", "Vanishing second derivative");
44 FlagDefinition const SdssCentroidAlgorithm::ALMOST_NO_SECOND_DERIVATIVE = flagDefinitions.add("flag_almostNoSecondDerivative", "Almost vanishing second derivative");
45 FlagDefinition const SdssCentroidAlgorithm::NOT_AT_MAXIMUM = flagDefinitions.add("flag_notAtMaximum", "Object is not at a maximum");
46 
48  return flagDefinitions;
49 }
50 
51 
52 namespace {
53 
54 /************************************************************************************************************/
55 
56 float const AMPAST4 = 1.33; // amplitude of `4th order' corr compared to theory
57 
58 /*
59  * Do the Gaussian quartic interpolation for the position
60  * of the maximum for three equally spaced values vm,v0,vp, assumed to be at
61  * abscissae -1,0,1; the answer is returned as *cen
62  *
63  * Return 0 is all is well, otherwise 1
64  */
65 static int inter4(float vm, float v0, float vp, float *cen, bool negative=false) {
66  float const sp = v0 - vp;
67  float const sm = v0 - vm;
68  float const d2 = sp + sm;
69  float const s = 0.5*(vp - vm);
70 
71  if ((!negative && (d2 <= 0.0f || v0 <= 0.0f)) ||
72  ( negative && (d2 >= 0.0f || v0 >= 0.0f))) {
73  return(1);
74  }
75 
76  *cen = s/d2*(1.0 + AMPAST4*sp*sm/(d2*v0));
77 
78  return fabs(*cen) < 1 ? 0 : 1;
79 }
80 
81 /*****************************************************************************/
82 /*
83  * Calculate error in centroid
84  */
85 float astrom_errors(float skyVar, // variance of pixels at the sky level
86  float sourceVar, // variance in peak due to excess counts over sky
87  float A, // abs(peak value in raw image)
88  float tau2, // Object is N(0,tau2)
89  float As, // abs(peak value in smoothed image)
90  float s, // slope across central pixel
91  float d, // curvature at central pixel
92  float sigma, // width of smoothing filter
93  int quarticBad) { // was quartic estimate bad?
94 
95  float const k = quarticBad ? 0 : AMPAST4; /* quartic correction coeff */
96  float const sigma2 = sigma*sigma; /* == sigma^2 */
97  float sVar, dVar; /* variances of s and d */
98  float xVar; /* variance of centroid, x */
99 
100  if (fabs(As) < std::numeric_limits<float>::min() ||
101  fabs(d) < std::numeric_limits<float>::min()) {
102  return(1e3);
103  }
104 
105  if (sigma <= 0) { /* no smoothing; no covariance */
106  sVar = 0.5*skyVar; /* due to sky */
107  dVar = 6*skyVar;
108 
109  sVar += 0.5*sourceVar*exp(-1/(2*tau2));
110  dVar += sourceVar*(4*exp(-1/(2*tau2)) + 2*exp(-1/(2*tau2)));
111  } else { /* smoothed */
112  sVar = skyVar/(8*lsst::afw::geom::PI*sigma2)*(1 - exp(-1/sigma2));
113  dVar = skyVar/(2*lsst::afw::geom::PI*sigma2)*(3 - 4*exp(-1/(4*sigma2)) + exp(-1/sigma2));
114 
115  sVar += sourceVar/(12*lsst::afw::geom::PI*sigma2)*(exp(-1/(3*sigma2)) - exp(-1/sigma2));
116  dVar += sourceVar/(3*lsst::afw::geom::PI*sigma2)*(2 - 3*exp(-1/(3*sigma2)) + exp(-1/sigma2));
117  }
118 
119  xVar = sVar*pow(1/d + k/(4*As)*(1 - 12*s*s/(d*d)), 2) +
120  dVar*pow(s/(d*d) - k/(4*As)*8*s*s/(d*d*d), 2);
121 
122  return(xVar >= 0 ? sqrt(xVar) : NAN);
123 }
124 
125 /************************************************************************************************************/
126 /*
127  * Estimate the position of an object, assuming we know that it's approximately the size of the PSF
128  */
129 
130 template<typename ImageXy_locatorT, typename VarImageXy_locatorT>
131 void doMeasureCentroidImpl(double *xCenter, // output; x-position of object
132  double *dxc, // output; error in xCenter
133  double *yCenter, // output; y-position of object
134  double *dyc, // output; error in yCenter
135  double *sizeX2, double *sizeY2, // output; object widths^2 in x and y directions
136  ImageXy_locatorT im, // Locator for the pixel values
137  VarImageXy_locatorT vim, // Locator for the image containing the variance
138  double smoothingSigma, // Gaussian sigma of already-applied smoothing filter
139  FlagHandler flagHandler
140  )
141 {
142  /*
143  * find a first quadratic estimate
144  */
145  double const d2x = 2*im(0, 0) - im(-1, 0) - im(1, 0);
146  double const d2y = 2*im(0, 0) - im( 0, -1) - im(0, 1);
147  double const sx = 0.5*(im(1, 0) - im(-1, 0));
148  double const sy = 0.5*(im(0, 1) - im( 0, -1));
149 
150  if (d2x == 0.0 || d2y == 0.0) {
151  throw LSST_EXCEPT(
152  MeasurementError,
155  );
156  }
157  if (d2x < 0.0 || d2y < 0.0) {
158  throw LSST_EXCEPT(
159  MeasurementError,
161  (boost::format(": d2I/dx2, d2I/dy2 = %g %g") % d2x % d2y).str(),
163  );
164  }
165 
166  double const dx0 = sx/d2x;
167  double const dy0 = sy/d2y; // first guess
168 
169  if (fabs(dx0) > 10.0 || fabs(dy0) > 10.0) {
170  throw LSST_EXCEPT(
171  MeasurementError,
173  (boost::format(": sx, d2x, sy, d2y = %f %f %f %f") % sx % d2x % sy % d2y).str(),
175  );
176  }
177 
178  double vpk = im(0, 0) + 0.5*(sx*dx0 + sy*dy0); // height of peak in image
179  if (vpk < 0) {
180  vpk = -vpk;
181  }
182 /*
183  * now evaluate maxima on stripes
184  */
185  float m0x = 0, m1x = 0, m2x = 0;
186  float m0y = 0, m1y = 0, m2y = 0;
187 
188  int quarticBad = 0;
189  quarticBad += inter4(im(-1, -1), im( 0, -1), im( 1, -1), &m0x);
190  quarticBad += inter4(im(-1, 0), im( 0, 0), im( 1, 0), &m1x);
191  quarticBad += inter4(im(-1, 1), im( 0, 1), im( 1, 1), &m2x);
192 
193  quarticBad += inter4(im(-1, -1), im(-1, 0), im(-1, 1), &m0y);
194  quarticBad += inter4(im( 0, -1), im( 0, 0), im( 0, 1), &m1y);
195  quarticBad += inter4(im( 1, -1), im( 1, 0), im( 1, 1), &m2y);
196 
197  double xc, yc; // position of maximum
198  double sigmaX2, sigmaY2; // widths^2 in x and y of smoothed object
199 
200  if (quarticBad) { // >= 1 quartic interpolator is bad
201  xc = dx0;
202  yc = dy0;
203  sigmaX2 = vpk/d2x; // widths^2 in x
204  sigmaY2 = vpk/d2y; // and y
205  } else {
206  double const smx = 0.5*(m2x - m0x);
207  double const smy = 0.5*(m2y - m0y);
208  double const dm2x = m1x - 0.5*(m0x + m2x);
209  double const dm2y = m1y - 0.5*(m0y + m2y);
210  double const dx = m1x + dy0*(smx - dy0*dm2x); // first quartic approx
211  double const dy = m1y + dx0*(smy - dx0*dm2y);
212  double const dx4 = m1x + dy*(smx - dy*dm2x); // second quartic approx
213  double const dy4 = m1y + dx*(smy - dx*dm2y);
214 
215  xc = dx4;
216  yc = dy4;
217  sigmaX2 = vpk/d2x - (1 + 6*dx0*dx0)/4; // widths^2 in x
218  sigmaY2 = vpk/d2y - (1 + 6*dy0*dy0)/4; // and y
219  }
220  /*
221  * Now for the errors.
222  */
223  float tauX2 = sigmaX2; // width^2 of _un_ smoothed object
224  float tauY2 = sigmaY2;
225  tauX2 -= smoothingSigma*smoothingSigma; // correct for smoothing
226  tauY2 -= smoothingSigma*smoothingSigma;
227 
228  if (tauX2 <= smoothingSigma*smoothingSigma) { // problem; sigmaX2 must be bad
229  tauX2 = smoothingSigma*smoothingSigma;
230  }
231  if (tauY2 <= smoothingSigma*smoothingSigma) { // sigmaY2 must be bad
232  tauY2 = smoothingSigma*smoothingSigma;
233  }
234 
235  float const skyVar = (vim(-1, -1) + vim( 0, -1) + vim( 1, -1) +
236  vim(-1, 0) + vim( 1, 0) +
237  vim(-1, 1) + vim( 0, 1) + vim( 1, 1))/8.0; // Variance in sky
238  float const sourceVar = vim(0, 0); // extra variance of peak due to its photons
239  float const A = vpk*sqrt((sigmaX2/tauX2)*(sigmaY2/tauY2)); // peak of Unsmoothed object
240 
241  *xCenter = xc;
242  *yCenter = yc;
243 
244  *dxc = astrom_errors(skyVar, sourceVar, A, tauX2, vpk, sx, d2x, fabs(smoothingSigma), quarticBad);
245  *dyc = astrom_errors(skyVar, sourceVar, A, tauY2, vpk, sy, d2y, fabs(smoothingSigma), quarticBad);
246 
247  *sizeX2 = tauX2; // return the estimates of the (object size)^2
248  *sizeY2 = tauY2;
249 }
250 
251 template<typename MaskedImageXy_locatorT>
252 void doMeasureCentroidImpl(double *xCenter, // output; x-position of object
253  double *dxc, // output; error in xCenter
254  double *yCenter, // output; y-position of object
255  double *dyc, // output; error in yCenter
256  double *sizeX2, double *sizeY2, // output; object widths^2 in x and y directions
257  double *peakVal, // output; peak of object
258  MaskedImageXy_locatorT mim, // Locator for the pixel values
259  double smoothingSigma, // Gaussian sigma of already-applied smoothing filter
260  bool negative,
261  FlagHandler flagHandler
262  )
263 {
264  /*
265  * find a first quadratic estimate
266  */
267  double const d2x = 2*mim.image(0, 0) - mim.image(-1, 0) - mim.image(1, 0);
268  double const d2y = 2*mim.image(0, 0) - mim.image( 0, -1) - mim.image(0, 1);
269  double const sx = 0.5*(mim.image(1, 0) - mim.image(-1, 0));
270  double const sy = 0.5*(mim.image(0, 1) - mim.image( 0, -1));
271 
272  if (d2x == 0.0 || d2y == 0.0) {
273  throw LSST_EXCEPT(
274  MeasurementError,
277  );
278  }
279  if ((!negative && (d2x < 0.0 || d2y < 0.0)) ||
280  ( negative && (d2x > 0.0 || d2y > 0.0))) {
281  throw LSST_EXCEPT(
282  MeasurementError,
284  (boost::format(": d2I/dx2, d2I/dy2 = %g %g") % d2x % d2y).str(),
286  );
287  }
288 
289  double const dx0 = sx/d2x;
290  double const dy0 = sy/d2y; // first guess
291 
292  if (fabs(dx0) > 10.0 || fabs(dy0) > 10.0) {
293  throw LSST_EXCEPT(
294  MeasurementError,
296  (boost::format(": sx, d2x, sy, d2y = %f %f %f %f") % sx % d2x % sy % d2y).str(),
298  );
299  }
300 
301  double vpk = mim.image(0, 0) + 0.5*(sx*dx0 + sy*dy0); // height of peak in image
302  //if (vpk < 0) {
303  // vpk = -vpk;
304  //}
305 /*
306  * now evaluate maxima on stripes
307  */
308  float m0x = 0, m1x = 0, m2x = 0;
309  float m0y = 0, m1y = 0, m2y = 0;
310 
311  int quarticBad = 0;
312  quarticBad += inter4(mim.image(-1, -1), mim.image( 0, -1), mim.image( 1, -1), &m0x, negative);
313  quarticBad += inter4(mim.image(-1, 0), mim.image( 0, 0), mim.image( 1, 0), &m1x, negative);
314  quarticBad += inter4(mim.image(-1, 1), mim.image( 0, 1), mim.image( 1, 1), &m2x, negative);
315 
316  quarticBad += inter4(mim.image(-1, -1), mim.image(-1, 0), mim.image(-1, 1), &m0y, negative);
317  quarticBad += inter4(mim.image( 0, -1), mim.image( 0, 0), mim.image( 0, 1), &m1y, negative);
318  quarticBad += inter4(mim.image( 1, -1), mim.image( 1, 0), mim.image( 1, 1), &m2y, negative);
319 
320  double xc, yc; // position of maximum
321  double sigmaX2, sigmaY2; // widths^2 in x and y of smoothed object
322 
323  if (quarticBad) { // >= 1 quartic interpolator is bad
324  xc = dx0;
325  yc = dy0;
326  sigmaX2 = vpk/d2x; // widths^2 in x
327  sigmaY2 = vpk/d2y; // and y
328  } else {
329  double const smx = 0.5*(m2x - m0x);
330  double const smy = 0.5*(m2y - m0y);
331  double const dm2x = m1x - 0.5*(m0x + m2x);
332  double const dm2y = m1y - 0.5*(m0y + m2y);
333  double const dx = m1x + dy0*(smx - dy0*dm2x); // first quartic approx
334  double const dy = m1y + dx0*(smy - dx0*dm2y);
335  double const dx4 = m1x + dy*(smx - dy*dm2x); // second quartic approx
336  double const dy4 = m1y + dx*(smy - dx*dm2y);
337 
338  xc = dx4;
339  yc = dy4;
340  sigmaX2 = vpk/d2x - (1 + 6*dx0*dx0)/4; // widths^2 in x
341  sigmaY2 = vpk/d2y - (1 + 6*dy0*dy0)/4; // and y
342  }
343  /*
344  * Now for the errors.
345  */
346  float tauX2 = sigmaX2; // width^2 of _un_ smoothed object
347  float tauY2 = sigmaY2;
348  tauX2 -= smoothingSigma*smoothingSigma; // correct for smoothing
349  tauY2 -= smoothingSigma*smoothingSigma;
350 
351  if (tauX2 <= smoothingSigma*smoothingSigma) { // problem; sigmaX2 must be bad
352  tauX2 = smoothingSigma*smoothingSigma;
353  }
354  if (tauY2 <= smoothingSigma*smoothingSigma) { // sigmaY2 must be bad
355  tauY2 = smoothingSigma*smoothingSigma;
356  }
357 
358  float const skyVar = (mim.variance(-1, -1) + mim.variance( 0, -1) + mim.variance( 1, -1) +
359  mim.variance(-1, 0) + mim.variance( 1, 0) +
360  mim.variance(-1, 1) + mim.variance( 0, 1) + mim.variance( 1, 1)
361  )/8.0; // Variance in sky
362  float const sourceVar = mim.variance(0, 0); // extra variance of peak due to its photons
363  float const A = vpk*sqrt((sigmaX2/tauX2)*(sigmaY2/tauY2)); // peak of Unsmoothed object
364 
365  *xCenter = xc;
366  *yCenter = yc;
367 
368  *dxc = astrom_errors(skyVar, sourceVar, A, tauX2, vpk, sx, d2x, fabs(smoothingSigma), quarticBad);
369  *dyc = astrom_errors(skyVar, sourceVar, A, tauY2, vpk, sy, d2y, fabs(smoothingSigma), quarticBad);
370 
371  *sizeX2 = tauX2; // return the estimates of the (object size)^2
372  *sizeY2 = tauY2;
373 
374  *peakVal = vpk;
375 }
376 
377 template<typename MaskedImageT>
378 std::pair<MaskedImageT, double>
379 smoothAndBinImage(CONST_PTR(lsst::afw::detection::Psf) psf,
380  int const x, const int y,
381  MaskedImageT const& mimage,
382  int binX, int binY,
383  FlagHandler _flagHandler)
384 {
385  lsst::afw::geom::Point2D const center(x + mimage.getX0(), y + mimage.getY0());
386  lsst::afw::geom::ellipses::Quadrupole const& shape = psf->computeShape(center);
387  double const smoothingSigma = shape.getDeterminantRadius();
388 #if 0
389  double const nEffective = psf->computeEffectiveArea(); // not implemented yet (#2821)
390 #else
391  double const nEffective = 4*M_PI*smoothingSigma*smoothingSigma; // correct for a Gaussian
392 #endif
393 
394  std::shared_ptr<lsst::afw::math::Kernel const> kernel = psf->getLocalKernel(center);
395  int const kWidth = kernel->getWidth();
396  int const kHeight = kernel->getHeight();
397 
398  lsst::afw::geom::BoxI bbox(lsst::afw::geom::Point2I(x - binX*(2 + kWidth/2), y - binY*(2 + kHeight/2)),
399  lsst::afw::geom::ExtentI(binX*(3 + kWidth + 1), binY*(3 + kHeight + 1)));
400 
401  // image to smooth, a shallow copy
402  PTR(MaskedImageT) subImage;
403  try {
404  subImage.reset(new MaskedImageT(mimage, bbox, lsst::afw::image::LOCAL));
405  } catch (pex::exceptions::LengthError & err) {
406  throw LSST_EXCEPT(
407  MeasurementError,
410  );
411  }
412  PTR(MaskedImageT) binnedImage = lsst::afw::math::binImage(*subImage, binX, binY, lsst::afw::math::MEAN);
413  binnedImage->setXY0(subImage->getXY0());
414  // image to smooth into, a deep copy.
415  MaskedImageT smoothedImage = MaskedImageT(*binnedImage, true);
416  assert(smoothedImage.getWidth()/2 == kWidth/2 + 2); // assumed by the code that uses smoothedImage
417  assert(smoothedImage.getHeight()/2 == kHeight/2 + 2);
418 
419  lsst::afw::math::convolve(smoothedImage, *binnedImage, *kernel, lsst::afw::math::ConvolutionControl());
420  *smoothedImage.getVariance() *= binX*binY*nEffective; // We want the per-pixel variance, so undo the
421  // effects of binning and smoothing
422 
423  return std::make_pair(smoothedImage, smoothingSigma);
424 }
425 
426 
427 } // end anonymous namespace
428 
430  Control const & ctrl,
431  std::string const & name,
432  afw::table::Schema & schema
433 ) : _ctrl(ctrl),
434  _centroidKey(
435  CentroidResultKey::addFields(schema, name, "centroid from Sdss Centroid algorithm", SIGMA_ONLY)
436  ),
437  _flagHandler(FlagHandler::addFields(schema, name,
438  getFlagDefinitions())),
439  _centroidExtractor(schema, name, true),
440  _centroidChecker(schema, name, ctrl.doFootprintCheck, ctrl.maxDistToPeak)
441 {
442 }
444  afw::table::SourceRecord & measRecord,
445  afw::image::Exposure<float> const & exposure
446 ) const {
447 
448  // get our current best guess about the centroid: either a centroider measurement or peak.
449  afw::geom::Point2D center = _centroidExtractor(measRecord, _flagHandler);
450  CentroidResult result;
451  result.x = center.getX();
452  result.y = center.getY();
453  measRecord.set(_centroidKey, result); // better than NaN
454 
455  typedef afw::image::Exposure<float>::MaskedImageT MaskedImageT;
456  typedef MaskedImageT::Image ImageT;
457  typedef MaskedImageT::Variance VarianceT;
458  bool negative = false;
459  try {
460  negative = measRecord.get(measRecord.getSchema().find<afw::table::Flag>("flags_negative").key);
461  } catch(pexExcept::Exception &e) {}
462 
463  MaskedImageT const& mimage = exposure.getMaskedImage();
464  ImageT const& image = *mimage.getImage();
465  CONST_PTR(lsst::afw::detection::Psf) psf = exposure.getPsf();
466 
467  int const x = image.positionToIndex(center.getX(), lsst::afw::image::X).first;
468  int const y = image.positionToIndex(center.getY(), lsst::afw::image::Y).first;
469 
470  if (!image.getBBox().contains(lsst::afw::geom::Extent2I(x,y) + image.getXY0())) {
471  throw LSST_EXCEPT(
473  EDGE.doc,
474  EDGE.number
475  );
476  }
477 
478  // Algorithm uses a least-squares fit (implemented via a convolution) to a symmetrized PSF model.
479  // If you don't have a Psf, you need to use another centroider, such as GaussianCentroider.
480  if (!psf) {
481  throw LSST_EXCEPT(
483  "SdssCentroid algorithm requires a Psf with every exposure"
484  );
485  }
486 
487  int binX = 1;
488  int binY = 1;
489  double xc=0., yc=0., dxc=0., dyc=0.; // estimated centre and error therein
490  for(int binsize = 1; binsize <= _ctrl.binmax; binsize *= 2) {
491  std::pair<MaskedImageT, double> result = smoothAndBinImage(psf, x, y, mimage, binX, binY, _flagHandler);
492  MaskedImageT const smoothedImage = result.first;
493  double const smoothingSigma = result.second;
494 
495  MaskedImageT::xy_locator mim = smoothedImage.xy_at(smoothedImage.getWidth()/2,
496  smoothedImage.getHeight()/2);
497 
498  double sizeX2, sizeY2; // object widths^2 in x and y directions
499  double peakVal; // peak intensity in image
500 
501  doMeasureCentroidImpl(&xc, &dxc, &yc, &dyc, &sizeX2, &sizeY2, &peakVal, mim,
502  smoothingSigma, negative, _flagHandler);
503 
504  if(binsize > 1) {
505  // dilate from the lower left corner of central pixel
506  xc = (xc + 0.5)*binX - 0.5;
507  dxc *= binX;
508  sizeX2 *= binX*binX;
509 
510  yc = (yc + 0.5)*binY - 0.5;
511  dyc *= binY;
512  sizeY2 *= binY*binY;
513  }
514 
515  xc += x; // xc, yc are measured relative to pixel (x, y)
516  yc += y;
517 
518  double const fac = _ctrl.wfac*(1 + smoothingSigma*smoothingSigma);
519  double const facX2 = fac*binX*binX;
520  double const facY2 = fac*binY*binY;
521 
522  if (sizeX2 < facX2 && ::pow(xc - x, 2) < facX2 &&
523  sizeY2 < facY2 && ::pow(yc - y, 2) < facY2) {
524  if (binsize > 1 || _ctrl.peakMin < 0.0 || peakVal > _ctrl.peakMin) {
525  break;
526  }
527  }
528 
529  if (sizeX2 >= facX2 || ::pow(xc - x, 2) >= facX2) {
530  binX *= 2;
531  }
532  if (sizeY2 >= facY2 || ::pow(yc - y, 2) >= facY2) {
533  binY *= 2;
534  }
535  }
536  result.x = lsst::afw::image::indexToPosition(xc + image.getX0());
537  result.y = lsst::afw::image::indexToPosition(yc + image.getY0());
538 
539  result.xSigma = sqrt(dxc*dxc);
540  result.ySigma = sqrt(dyc*dyc);
541  measRecord.set(_centroidKey, result);
542  _centroidChecker(measRecord);
543 }
544 
545 
546 void SdssCentroidAlgorithm::fail(afw::table::SourceRecord & measRecord, MeasurementError * error) const {
547  _flagHandler.handleFailure(measRecord, error);
548 }
549 
551  Control const & ctrl,
552  std::string const & name,
553  afw::table::SchemaMapper & mapper
554 ) :
555  CentroidTransform{name, mapper}
556 {
557  for (std::size_t i = 0; i < SdssCentroidAlgorithm::getFlagDefinitions().size(); i++) {
559  if (flag == SdssCentroidAlgorithm::FAILURE) continue;
560  if (mapper.getInputSchema().getNames().count(
561  mapper.getInputSchema().join(name, flag.name)) == 0) continue;
562  afw::table::Key<afw::table::Flag> key = mapper.getInputSchema().find<afw::table::Flag>(
563  name + "_" + flag.name).key;
564  mapper.addMapping(key);
565  }
566 }
567 
568 }}} // end namespace lsst::meas::base
569 
int x
The Sdss Centroid Algorithm.
Definition: SdssCentroid.h:68
ErrElement ySigma
1-Sigma uncertainty on y (sqrt of variance)
virtual void fail(afw::table::SourceRecord &measRecord, MeasurementError *error=nullptr) const
Handle an exception thrown by the current algorithm by setting flags in the given record...
Simple class used to define and document flags The name and doc constitute the identity of the FlagDe...
Definition: FlagHandler.h:38
Only the diagonal elements of the covariance matrix are provided.
Definition: constants.h:43
A reusable struct for centroid measurements.
ErrElement xSigma
1-Sigma uncertainty on x (sqrt of variance)
CentroidElement x
x (column) coordinate of the measured position
A C++ control class to handle SdssCentroidAlgorithm&#39;s configuration.
Definition: SdssCentroid.h:46
int binmax
&quot;maximum allowed binning&quot; ;
Definition: SdssCentroid.h:49
static FlagDefinition const ALMOST_NO_SECOND_DERIVATIVE
Definition: SdssCentroid.h:76
Exception to be thrown when a measurement algorithm experiences a known failure mode.
Definition: exceptions.h:48
static FlagDefinition const NOT_AT_MAXIMUM
Definition: SdssCentroid.h:77
static FlagDefinitionList const & getFlagDefinitions()
Definition: SdssCentroid.cc:47
Exception to be thrown when a measurement algorithm experiences a fatal error.
Definition: exceptions.h:83
Utility class for handling flag fields that indicate the failure modes of an algorithm.
Definition: FlagHandler.h:156
_centroidExtractor(schema, name)
Definition: ApertureFlux.cc:85
static FlagDefinition const EDGE
Definition: SdssCentroid.h:74
static FlagDefinition const NO_SECOND_DERIVATIVE
Definition: SdssCentroid.h:75
double peakMin
&quot;if the peak&#39;s less than this insist on binning at least once&quot; ;
Definition: SdssCentroid.h:50
double wfac
&quot;fiddle factor for adjusting the binning&quot; ;
Definition: SdssCentroid.h:51
virtual void measure(afw::table::SourceRecord &measRecord, afw::image::Exposure< float > const &exposure) const
Called to measure a single child source in an image.
A FunctorKey for CentroidResult.
CentroidElement y
y (row) coordinate of the measured position
afw::geom::Point< CentroidElement, 2 > Centroid
Definition: constants.h:56
This implements the SdssCentroid algorithm within the meas_base measurement framework.
vector-type utility class to build a collection of FlagDefinitions
Definition: FlagHandler.h:63
std::size_t size() const
return the current size (number of defined elements) of the collection
Definition: FlagHandler.h:133
static FlagDefinition const FAILURE
Definition: SdssCentroid.h:73
Base for centroid measurement transformations.
SdssCentroidTransform(Control const &ctrl, std::string const &name, afw::table::SchemaMapper &mapper)
int y