lsst.meas.base  14.0-14-gdfadacb+2
SdssShape.cc
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 /*
3  * LSST Data Management System
4  * Copyright 2008-2016 AURA/LSST.
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 <https://www.lsstcorp.org/LegalNotices/>.
22  */
23 
24 #include <cmath>
25 #include <tuple>
26 
27 #include "boost/tuple/tuple.hpp"
28 #include "Eigen/LU"
29 #include "lsst/afw/image.h"
30 #include "lsst/afw/detection/Psf.h"
31 #include "lsst/afw/geom/Angle.h"
32 #include "lsst/afw/geom/ellipses.h"
33 #include "lsst/afw/table/Source.h"
34 
37 
38 namespace pexPolicy = lsst::pex::policy;
40 namespace afwDet = lsst::afw::detection;
41 namespace afwImage = lsst::afw::image;
42 namespace afwGeom = lsst::afw::geom;
43 namespace afwTable = lsst::afw::table;
44 
45 namespace lsst { namespace meas { namespace base {
46 namespace {
47 FlagDefinitionList flagDefinitions;
48 } // end anonymous
49 
50 FlagDefinition const SdssShapeAlgorithm::FAILURE = flagDefinitions.addFailureFlag();
51 FlagDefinition const SdssShapeAlgorithm::UNWEIGHTED_BAD = flagDefinitions.add("flag_unweightedBad", "Both weighted and unweighted moments were invalid");
52 FlagDefinition const SdssShapeAlgorithm::UNWEIGHTED = flagDefinitions.add("flag_unweighted", "Weighted moments converged to an invalid value; using unweighted moments");
53 FlagDefinition const SdssShapeAlgorithm::SHIFT = flagDefinitions.add("flag_shift", "centroid shifted by more than the maximum allowed amount");
54 FlagDefinition const SdssShapeAlgorithm::MAXITER = flagDefinitions.add("flag_maxIter", "Too many iterations in adaptive moments");
55 FlagDefinition const SdssShapeAlgorithm::PSF_SHAPE_BAD = flagDefinitions.add("flag_psf", "Failure in measuring PSF model shape");
56 
58  return flagDefinitions;
59 }
60 
61 
62 namespace { // anonymous
63 
64 typedef Eigen::Matrix<double,4,4,Eigen::DontAlign> Matrix4d;
65 
66 // Return multiplier that transforms I0 to a total flux
67 double computeFluxScale(SdssShapeResult const & result) {
68  /*
69  * The shape is an ellipse that's axis-aligned in (u, v) [<uv> = 0] after rotation by theta:
70  * <x^2> + <y^2> = <u^2> + <v^2>
71  * <x^2> - <y^2> = cos(2 theta)*(<u^2> - <v^2>)
72  * 2*<xy> = sin(2 theta)*(<u^2> - <v^2>)
73  */
74  double const Mxx = result.xx; // <x^2>
75  double const Mxy = result.xy; // <xy>
76  double const Myy = result.yy; // <y^2>
77 
78  double const Muu_p_Mvv = Mxx + Myy; // <u^2> + <v^2>
79  double const Muu_m_Mvv = ::sqrt(::pow(Mxx - Myy, 2) + 4*::pow(Mxy, 2)); // <u^2> - <v^2>
80  double const Muu = 0.5*(Muu_p_Mvv + Muu_m_Mvv);
81  double const Mvv = 0.5*(Muu_p_Mvv - Muu_m_Mvv);
82 
83  return lsst::afw::geom::TWOPI * ::sqrt(Muu*Mvv);
84 }
85 
86 /*****************************************************************************/
87 /*
88  * Error analysis, courtesy of David Johnston, University of Chicago
89  */
90 /*
91  * This function takes the 4 Gaussian parameters A, sigmaXXW and the
92  * sky variance and fills in the Fisher matrix from the least squares fit.
93  *
94  * Following "Numerical Recipes in C" section 15.5, it ignores the 2nd
95  * derivative parts and so the fisher matrix is just a function of these
96  * best fit model parameters. The components are calculated analytically.
97  */
98 Matrix4d
99 calc_fisher(SdssShapeResult const& shape, // the Shape that we want the the Fisher matrix for
100  float bkgd_var // background variance level for object
101 ) {
102  float const A = shape.flux; // amplitude; will be converted to flux later
103  float const sigma11W = shape.xx;
104  float const sigma12W = shape.xy;
105  float const sigma22W = shape.yy;
106 
107  double const D = sigma11W*sigma22W - sigma12W*sigma12W;
108 
111  "Determinant is too small calculating Fisher matrix");
112  }
113 /*
114  * a normalization factor
115  */
116  if (bkgd_var <= 0.0) {
118  (boost::format("Background variance must be positive (saw %g)") % bkgd_var).str());
119  }
120  double const F = afwGeom::PI*sqrt(D)/bkgd_var;
121 /*
122  * Calculate the 10 independent elements of the 4x4 Fisher matrix
123  */
124  Matrix4d fisher;
125 
126  double fac = F*A/(4.0*D);
127  fisher(0, 0) = F;
128  fisher(0, 1) = fac*sigma22W;
129  fisher(1, 0) = fisher(0, 1);
130  fisher(0, 2) = fac*sigma11W;
131  fisher(2, 0) = fisher(0, 2);
132  fisher(0, 3) = -fac*2*sigma12W;
133  fisher(3, 0) = fisher(0, 3);
134 
135  fac = 3.0*F*A*A/(16.0*D*D);
136  fisher(1, 1) = fac*sigma22W*sigma22W;
137  fisher(2, 2) = fac*sigma11W*sigma11W;
138  fisher(3, 3) = fac*4.0*(sigma12W*sigma12W + D/3.0);
139 
140  fisher(1, 2) = fisher(3, 3)/4.0;
141  fisher(2, 1) = fisher(1, 2);
142  fisher(1, 3) = fac*(-2*sigma22W*sigma12W);
143  fisher(3, 1) = fisher(1, 3);
144  fisher(2, 3) = fac*(-2*sigma11W*sigma12W);
145  fisher(3, 2) = fisher(2, 3);
146 
147  return fisher;
148 }
149 //
150 // Here's a class to allow us to get the Image and variance from an Image or MaskedImage
151 //
152 template<typename ImageT> // general case
153 struct ImageAdaptor {
154  typedef ImageT Image;
155 
156  static bool const hasVariance = false;
157 
158  Image const& getImage(ImageT const& image) const {
159  return image;
160  }
161 
162  double getVariance(ImageT const&, int, int) {
164  }
165 };
166 
167 template<typename T> // specialise to a MaskedImage
168 struct ImageAdaptor<afwImage::MaskedImage<T> > {
169  typedef typename afwImage::MaskedImage<T>::Image Image;
170 
171  static bool const hasVariance = true;
172 
173  Image const& getImage(afwImage::MaskedImage<T> const& mimage) const {
174  return *mimage.getImage();
175  }
176 
177  double getVariance(afwImage::MaskedImage<T> const& mimage, int ix, int iy) {
178  return mimage.at(ix, iy).variance();
179  }
180 };
181 
183 std::tuple<std::pair<bool, double>, double, double, double>
184 getWeights(double sigma11, double sigma12, double sigma22) {
185  double const NaN = std::numeric_limits<double>::quiet_NaN();
186  if (std::isnan(sigma11) || std::isnan(sigma12) || std::isnan(sigma22)) {
187  return std::make_tuple(std::make_pair(false, NaN), NaN, NaN, NaN);
188  }
189  double const det = sigma11*sigma22 - sigma12*sigma12; // determinant of sigmaXX matrix
190  if (std::isnan(det) || det < std::numeric_limits<float>::epsilon()) { // a suitably small number
191  return std::make_tuple(std::make_pair(false, det), NaN, NaN, NaN);
192  }
193  return std::make_tuple(std::make_pair(true, det), sigma22/det, -sigma12/det, sigma11/det);
194 }
195 
197 bool shouldInterp(double sigma11, double sigma22, double det) {
198  float const xinterp = 0.25; // I.e. 0.5*0.5
199  return (sigma11 < xinterp || sigma22 < xinterp || det < xinterp*xinterp);
200 }
201 
202 // Decide on the bounding box for the region to examine while calculating the adaptive moments
203 // This routine will work in either LOCAL or PARENT coordinates (but of course which you pass
204 // determine which you will get back).
205 afw::geom::Box2I computeAdaptiveMomentsBBox(
206  afw::geom::Box2I const & bbox, // full image bbox
207  afw::geom::Point2D const & center, // centre of object
208  double sigma11_w, // quadratic moments of the
209  double , // weighting function
210  double sigma22_w, // xx, xy, and yy
211  double maxRadius = 1000 // Maximum radius of area to use
212 ) {
213  double radius = std::min(4*std::sqrt(std::max(sigma11_w, sigma22_w)), maxRadius);
214  afw::geom::Extent2D offset(radius, radius);
215  afw::geom::Box2I result(afw::geom::Box2D(center - offset, center + offset));
216  result.clip(bbox);
217  return result;
218 }
219 
220 /*****************************************************************************/
221 /*
222  * Calculate weighted moments of an object up to 2nd order
223  */
224 template<bool fluxOnly, typename ImageT>
225 static int
226 calcmom(ImageT const& image, // the image data
227  float xcen, float ycen, // centre of object
228  lsst::afw::geom::BoxI bbox, // bounding box to consider
229  float bkgd, // data's background level
230  bool interpflag, // interpolate within pixels?
231  double w11, double w12, double w22, // weights
232  double *pI0, // amplitude of fit
233  double *psum, // sum w*I (if !NULL)
234  double *psumx, double *psumy, // sum [xy]*w*I (if !fluxOnly)
235  double *psumxx, double *psumxy, double *psumyy, // sum [xy]^2*w*I (if !fluxOnly)
236  double *psums4, // sum w*I*weight^2 (if !fluxOnly && !NULL)
237  bool negative = false
238  )
239 {
240 
241  float tmod, ymod;
242  float X, Y; // sub-pixel interpolated [xy]
243  float weight;
244  float tmp;
245  double sum, sumx, sumy, sumxx, sumyy, sumxy, sums4;
246 #define RECALC_W 0 // estimate sigmaXX_w within BBox?
247 #if RECALC_W
248  double wsum, wsumxx, wsumxy, wsumyy;
249 
250  wsum = wsumxx = wsumxy = wsumyy = 0;
251 #endif
252 
253  assert(w11 >= 0); // i.e. it was set
254  if (fabs(w11) > 1e6 || fabs(w12) > 1e6 || fabs(w22) > 1e6) {
255  return(-1);
256  }
257 
258  sum = sumx = sumy = sumxx = sumxy = sumyy = sums4 = 0;
259 
260  int const ix0 = bbox.getMinX(); // corners of the box being analyzed
261  int const ix1 = bbox.getMaxX();
262  int const iy0 = bbox.getMinY(); // corners of the box being analyzed
263  int const iy1 = bbox.getMaxY();
264 
265  if (ix0 < 0 || ix1 >= image.getWidth() || iy0 < 0 || iy1 >= image.getHeight()) {
266  return -1;
267  }
268 
269  for (int i = iy0; i <= iy1; ++i) {
270  typename ImageT::x_iterator ptr = image.x_at(ix0, i);
271  float const y = i - ycen;
272  float const y2 = y*y;
273  float const yl = y - 0.375;
274  float const yh = y + 0.375;
275  for (int j = ix0; j <= ix1; ++j, ++ptr) {
276  float x = j - xcen;
277  if (interpflag) {
278  float const xl = x - 0.375;
279  float const xh = x + 0.375;
280 
281  float expon = xl*xl*w11 + yl*yl*w22 + 2.0*xl*yl*w12;
282  tmp = xh*xh*w11 + yh*yh*w22 + 2.0*xh*yh*w12;
283  expon = (expon > tmp) ? expon : tmp;
284  tmp = xl*xl*w11 + yh*yh*w22 + 2.0*xl*yh*w12;
285  expon = (expon > tmp) ? expon : tmp;
286  tmp = xh*xh*w11 + yl*yl*w22 + 2.0*xh*yl*w12;
287  expon = (expon > tmp) ? expon : tmp;
288 
289  if (expon <= 9.0) {
290  tmod = *ptr - bkgd;
291  for (Y = yl; Y <= yh; Y += 0.25) {
292  double const interpY2 = Y*Y;
293  for (X = xl; X <= xh; X += 0.25) {
294  double const interpX2 = X*X;
295  double const interpXy = X*Y;
296  expon = interpX2*w11 + 2*interpXy*w12 + interpY2*w22;
297  weight = std::exp(-0.5*expon);
298 
299  ymod = tmod*weight;
300  sum += ymod;
301  if (!fluxOnly) {
302  sumx += ymod*(X + xcen);
303  sumy += ymod*(Y + ycen);
304 #if RECALC_W
305  wsum += weight;
306 
307  tmp = interpX2*weight;
308  wsumxx += tmp;
309  sumxx += tmod*tmp;
310 
311  tmp = interpXy*weight;
312  wsumxy += tmp;
313  sumxy += tmod*tmp;
314 
315  tmp = interpY2*weight;
316  wsumyy += tmp;
317  sumyy += tmod*tmp;
318 #else
319  sumxx += interpX2*ymod;
320  sumxy += interpXy*ymod;
321  sumyy += interpY2*ymod;
322 #endif
323  sums4 += expon*expon*ymod;
324  }
325  }
326  }
327  }
328  } else {
329  float x2 = x*x;
330  float xy = x*y;
331  float expon = x2*w11 + 2*xy*w12 + y2*w22;
332 
333  if (expon <= 14.0) {
334  weight = std::exp(-0.5*expon);
335  tmod = *ptr - bkgd;
336  ymod = tmod*weight;
337  sum += ymod;
338  if (!fluxOnly) {
339  sumx += ymod*j;
340  sumy += ymod*i;
341 #if RECALC_W
342  wsum += weight;
343 
344  tmp = x2*weight;
345  wsumxx += tmp;
346  sumxx += tmod*tmp;
347 
348  tmp = xy*weight;
349  wsumxy += tmp;
350  sumxy += tmod*tmp;
351 
352  tmp = y2*weight;
353  wsumyy += tmp;
354  sumyy += tmod*tmp;
355 #else
356  sumxx += x2*ymod;
357  sumxy += xy*ymod;
358  sumyy += y2*ymod;
359 #endif
360  sums4 += expon*expon*ymod;
361  }
362  }
363  }
364  }
365  }
366 
367 
368  std::tuple<std::pair<bool, double>, double, double, double> const weights = getWeights(w11, w12, w22);
369  double const detW = std::get<1>(weights)*std::get<3>(weights) - std::pow(std::get<2>(weights), 2);
370  *pI0 = sum/(afwGeom::PI*sqrt(detW));
371  if (psum) {
372  *psum = sum;
373  }
374  if (!fluxOnly) {
375  *psumx = sumx;
376  *psumy = sumy;
377  *psumxx = sumxx;
378  *psumxy = sumxy;
379  *psumyy = sumyy;
380  if (psums4 != NULL) {
381  *psums4 = sums4;
382  }
383  }
384 
385 #if RECALC_W
386  if (wsum > 0 && !fluxOnly) {
387  double det = w11*w22 - w12*w12;
388  wsumxx /= wsum;
389  wsumxy /= wsum;
390  wsumyy /= wsum;
391  printf("%g %g %g %g %g %g\n", w22/det, -w12/det, w11/det, wsumxx, wsumxy, wsumyy);
392  }
393 #endif
394 
395  if (negative) {
396  return (fluxOnly || (sum < 0 && sumxx < 0 && sumyy < 0)) ? 0 : -1;
397  } else {
398  return (fluxOnly || (sum > 0 && sumxx > 0 && sumyy > 0)) ? 0 : -1;
399  }
400 }
401 
402 /*
403  * Workhorse for adaptive moments
404  *
405  * All inputs are expected to be in LOCAL image coordinates
406  */
407 template<typename ImageT>
408 bool getAdaptiveMoments(ImageT const& mimage, double bkgd, double xcen, double ycen, double shiftmax,
409  SdssShapeResult *shape, int maxIter, float tol1, float tol2, bool negative)
410 {
411  double I0 = 0; // amplitude of best-fit Gaussian
412  double sum; // sum of intensity*weight
413  double sumx, sumy; // sum ((int)[xy])*intensity*weight
414  double sumxx, sumxy, sumyy; // sum {x^2,xy,y^2}*intensity*weight
415  double sums4; // sum intensity*weight*exponent^2
416  float const xcen0 = xcen; // initial centre
417  float const ycen0 = ycen; // of object
418 
419  double sigma11W = 1.5; // quadratic moments of the
420  double sigma12W = 0.0; // weighting fcn;
421  double sigma22W = 1.5; // xx, xy, and yy
422 
423  double w11 = -1, w12 = -1, w22 = -1; // current weights for moments; always set when iter == 0
424  float e1_old = 1e6, e2_old = 1e6; // old values of shape parameters e1 and e2
425  float sigma11_ow_old = 1e6; // previous version of sigma11_ow
426 
427  typename ImageAdaptor<ImageT>::Image const &image = ImageAdaptor<ImageT>().getImage(mimage);
428 
429  if (std::isnan(xcen) || std::isnan(ycen)) {
430  // Can't do anything
431  shape->flags[SdssShapeAlgorithm::UNWEIGHTED_BAD.number] = true;
432  return false;
433  }
434 
435  bool interpflag = false; // interpolate finer than a pixel?
437  int iter = 0; // iteration number
438  for (; iter < maxIter; iter++) {
439  bbox = computeAdaptiveMomentsBBox(image.getBBox(afw::image::LOCAL), afw::geom::Point2D(xcen, ycen),
440  sigma11W, sigma12W, sigma22W);
441  std::tuple<std::pair<bool, double>, double, double, double> weights =
442  getWeights(sigma11W, sigma12W, sigma22W);
443  if (!std::get<0>(weights).first) {
444  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
445  break;
446  }
447 
448  double const detW = std::get<0>(weights).second;
449 
450 #if 0 // this form was numerically unstable on my G4 powerbook
451  assert(detW >= 0.0);
452 #else
453  assert(sigma11W*sigma22W >= sigma12W*sigma12W - std::numeric_limits<float>::epsilon());
454 #endif
455 
456  {
457  const double ow11 = w11; // old
458  const double ow12 = w12; // values
459  const double ow22 = w22; // of w11, w12, w22
460 
461  w11 = std::get<1>(weights);
462  w12 = std::get<2>(weights);
463  w22 = std::get<3>(weights);
464 
465  if (shouldInterp(sigma11W, sigma22W, detW)) {
466  if (!interpflag) {
467  interpflag = true; // N.b.: stays set for this object
468  if (iter > 0) {
469  sigma11_ow_old = 1.e6; // force at least one more iteration
470  w11 = ow11;
471  w12 = ow12;
472  w22 = ow22;
473  iter--; // we didn't update wXX
474  }
475  }
476  }
477  }
478 
479  if (calcmom<false>(image, xcen, ycen, bbox, bkgd, interpflag, w11, w12, w22,
480  &I0, &sum, &sumx, &sumy, &sumxx, &sumxy, &sumyy, &sums4, negative) < 0) {
481  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
482  break;
483  }
484 
485 #if 0
486 /*
487  * Find new centre
488  *
489  * This is only needed if we update the centre; if we use the input position we've already done the work
490  */
491  xcen = sumx/sum;
492  ycen = sumy/sum;
493 #endif
494  shape->x = sumx/sum; // update centroid. N.b. we're not setting errors here
495  shape->y = sumy/sum;
496 
497  if (fabs(shape->x - xcen0) > shiftmax || fabs(shape->y - ycen0) > shiftmax) {
498  shape->flags[SdssShapeAlgorithm::SHIFT.number] = true;
499  }
500 /*
501  * OK, we have the centre. Proceed to find the second moments.
502  */
503  float const sigma11_ow = sumxx/sum; // quadratic moments of
504  float const sigma22_ow = sumyy/sum; // weight*object
505  float const sigma12_ow = sumxy/sum; // xx, xy, and yy
506 
507  if (sigma11_ow <= 0 || sigma22_ow <= 0) {
508  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
509  break;
510  }
511 
512  float const d = sigma11_ow + sigma22_ow; // current values of shape parameters
513  float const e1 = (sigma11_ow - sigma22_ow)/d;
514  float const e2 = 2.0*sigma12_ow/d;
515 /*
516  * Did we converge?
517  */
518  if (iter > 0 &&
519  fabs(e1 - e1_old) < tol1 && fabs(e2 - e2_old) < tol1 &&
520  fabs(sigma11_ow/sigma11_ow_old - 1.0) < tol2 ) {
521  break; // yes; we converged
522  }
523 
524  e1_old = e1;
525  e2_old = e2;
526  sigma11_ow_old = sigma11_ow;
527 /*
528  * Didn't converge, calculate new values for weighting function
529  *
530  * The product of two Gaussians is a Gaussian:
531  * <x^2 exp(-a x^2 - 2bxy - cy^2) exp(-Ax^2 - 2Bxy - Cy^2)> =
532  * <x^2 exp(-(a + A) x^2 - 2(b + B)xy - (c + C)y^2)>
533  * i.e. the inverses of the covariances matrices add.
534  *
535  * We know sigmaXX_ow and sigmaXXW, the covariances of the weighted object
536  * and of the weights themselves. We can estimate the object's covariance as
537  * sigmaXX_ow^-1 - sigmaXXW^-1
538  * and, as we want to find a set of weights with the _same_ covariance as the
539  * object we take this to be the an estimate of our correct weights.
540  *
541  * N.b. This assumes that the object is roughly Gaussian.
542  * Consider the object:
543  * O == delta(x + p) + delta(x - p)
544  * the covariance of the weighted object is equal to that of the unweighted
545  * object, and this prescription fails badly. If we detect this, we set
546  * the UNWEIGHTED flag, and calculate the UNweighted moments
547  * instead.
548  */
549  {
550  float n11, n12, n22; // elements of inverse of next guess at weighting function
551  float ow11, ow12, ow22; // elements of inverse of sigmaXX_ow
552 
553  std::tuple<std::pair<bool, double>, double, double, double> weights =
554  getWeights(sigma11_ow, sigma12_ow, sigma22_ow);
555  if (!std::get<0>(weights).first) {
556  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
557  break;
558  }
559 
560  ow11 = std::get<1>(weights);
561  ow12 = std::get<2>(weights);
562  ow22 = std::get<3>(weights);
563 
564  n11 = ow11 - w11;
565  n12 = ow12 - w12;
566  n22 = ow22 - w22;
567 
568  weights = getWeights(n11, n12, n22);
569  if (!std::get<0>(weights).first) {
570  // product-of-Gaussians assumption failed
571  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
572  break;
573  }
574 
575  sigma11W = std::get<1>(weights);
576  sigma12W = std::get<2>(weights);
577  sigma22W = std::get<3>(weights);
578  }
579 
580  if (sigma11W <= 0 || sigma22W <= 0) {
581  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
582  break;
583  }
584  }
585 
586  if (iter == maxIter) {
587  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
588  shape->flags[SdssShapeAlgorithm::MAXITER.number] = true;
589  }
590 
591  if (sumxx + sumyy == 0.0) {
592  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = true;
593  }
594 /*
595  * Problems; try calculating the un-weighted moments
596  */
597  if (shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number]) {
598  w11 = w22 = w12 = 0;
599  if (calcmom<false>(image, xcen, ycen, bbox, bkgd, interpflag, w11, w12, w22,
600  &I0, &sum, &sumx, &sumy, &sumxx, &sumxy, &sumyy, NULL, negative) < 0 ||
601  (!negative && sum <= 0) || (negative && sum >= 0)) {
602  shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number] = false;
603  shape->flags[SdssShapeAlgorithm::UNWEIGHTED_BAD.number] = true;
604 
605  if (sum > 0) {
606  shape->xx = 1/12.0; // a single pixel
607  shape->xy = 0.0;
608  shape->yy = 1/12.0;
609  }
610 
611  return false;
612  }
613 
614  sigma11W = sumxx/sum; // estimate of object moments
615  sigma12W = sumxy/sum; // usually, object == weight
616  sigma22W = sumyy/sum; // at this point
617  }
618 
619  shape->flux = I0;
620  shape->xx = sigma11W;
621  shape->xy = sigma12W;
622  shape->yy = sigma22W;
623 
624  if (shape->xx + shape->yy != 0.0) {
625  int const ix = lsst::afw::image::positionToIndex(xcen);
626  int const iy = lsst::afw::image::positionToIndex(ycen);
627 
628  if (ix >= 0 && ix < mimage.getWidth() && iy >= 0 && iy < mimage.getHeight()) {
629  float const bkgd_var =
630  ImageAdaptor<ImageT>().getVariance(mimage, ix, iy); // XXX Overestimate as it includes object
631 
632  if (bkgd_var > 0.0) { // NaN is not > 0.0
633  if (!(shape->flags[SdssShapeAlgorithm::UNWEIGHTED.number])) {
634  Matrix4d fisher = calc_fisher(*shape, bkgd_var); // Fisher matrix
635  Matrix4d cov = fisher.inverse();
636  // convention in afw::geom::ellipses is to order moments (xx, yy, xy),
637  // but the older algorithmic code uses (xx, xy, yy) - the order of
638  // indices here is not a bug.
639  shape->fluxSigma = std::sqrt(cov(0, 0));
640  shape->xxSigma = std::sqrt(cov(1, 1));
641  shape->xySigma = std::sqrt(cov(2, 2));
642  shape->yySigma = std::sqrt(cov(3, 3));
643  shape->flux_xx_Cov = cov(0, 1);
644  shape->flux_xy_Cov = cov(0, 2);
645  shape->flux_yy_Cov = cov(0, 3);
646  shape->xx_yy_Cov = cov(1, 3);
647  shape->xx_xy_Cov = cov(1, 2);
648  shape->yy_xy_Cov = cov(2, 3);
649  }
650  }
651  }
652  }
653 
654  return true;
655 }
656 
657 } // anonymous
658 
659 
661  flux_xx_Cov(std::numeric_limits<ErrElement>::quiet_NaN()),
662  flux_yy_Cov(std::numeric_limits<ErrElement>::quiet_NaN()),
663  flux_xy_Cov(std::numeric_limits<ErrElement>::quiet_NaN())
664 {}
665 
666 
668  afw::table::Schema & schema,
669  std::string const & name,
670  bool doMeasurePsf
671 ) {
672 
674  r._shapeResult = ShapeResultKey::addFields(schema, name, "elliptical Gaussian adaptive moments",
675  SIGMA_ONLY);
676  r._centroidResult = CentroidResultKey::addFields(schema, name, "elliptical Gaussian adaptive moments",
678  r._fluxResult = FluxResultKey::addFields(schema, name, "elliptical Gaussian adaptive moments");
679 
680  // Only include PSF shape fields if doMeasurePsf = True
681  if (doMeasurePsf) {
682  r._includePsf = true;
683  r._psfShapeResult = afwTable::QuadrupoleKey::addFields(
684  schema, schema.join(name, "psf"), "adaptive moments of the PSF model at the object position");
685  } else {
686  r._includePsf = false;
687  }
688 
689  r._flux_xx_Cov = schema.addField<ErrElement>(
690  schema.join(name, "flux", "xx", "Cov"),
691  (boost::format("uncertainty covariance between %s and %s")
692  % schema.join(name, "flux") % schema.join(name, "xx")).str(),
693  "count*pixel^2"
694  );
695  r._flux_yy_Cov = schema.addField<ErrElement>(
696  schema.join(name, "flux", "yy", "Cov"),
697  (boost::format("uncertainty covariance between %s and %s")
698  % schema.join(name, "flux") % schema.join(name, "yy")).str(),
699  "count*pixel^2"
700  );
701  r._flux_xy_Cov = schema.addField<ErrElement>(
702  schema.join(name, "flux", "xy", "Cov"),
703  (boost::format("uncertainty covariance between %s and %s")
704  % schema.join(name, "flux") % schema.join(name, "xy")).str(),
705  "count*pixel^2"
706  );
707 
708  // Skip the psf flag if not recording the PSF shape.
709  if (r._includePsf) {
710  r._flagHandler = FlagHandler::addFields(schema, name, SdssShapeAlgorithm::getFlagDefinitions());
711  }
712  else {
713  r._flagHandler = FlagHandler::addFields(schema, name, SdssShapeAlgorithm::getFlagDefinitions(),
714  {SdssShapeAlgorithm::PSF_SHAPE_BAD});
715  }
716  return r;
717 }
718 
720  _shapeResult(s),
721  _centroidResult(s),
722  _fluxResult(s),
723  _flux_xx_Cov(s["flux"]["xx"]["Cov"]),
724  _flux_yy_Cov(s["flux"]["yy"]["Cov"]),
725  _flux_xy_Cov(s["flux"]["xy"]["Cov"])
726 {
727  // The input SubSchema may optionally provide for a PSF.
728  try {
729  _psfShapeResult = afwTable::QuadrupoleKey(s["psf"]);
731  _includePsf = true;
732  } catch (pex::exceptions::NotFoundError& e) {
733  _flagHandler = FlagHandler(s, SdssShapeAlgorithm::getFlagDefinitions(), {SdssShapeAlgorithm::PSF_SHAPE_BAD});
734  _includePsf = false;
735  }
736 }
737 
739  SdssShapeResult result;
740  static_cast<ShapeResult&>(result) = record.get(_shapeResult);
741  static_cast<CentroidResult&>(result) = record.get(_centroidResult);
742  static_cast<FluxResult&>(result) = record.get(_fluxResult);
743  result.flux_xx_Cov = record.get(_flux_xx_Cov);
744  result.flux_yy_Cov = record.get(_flux_yy_Cov);
745  result.flux_xy_Cov = record.get(_flux_xy_Cov);
746  for (size_t n = 0; n < SdssShapeAlgorithm::N_FLAGS; ++n) {
747  if (n == SdssShapeAlgorithm::PSF_SHAPE_BAD.number && !_includePsf) continue;
748  result.flags[n] = _flagHandler.getValue(record, n);
749  }
750  return result;
751 }
752 
753 afwGeom::ellipses::Quadrupole SdssShapeResultKey::getPsfShape(afw::table::BaseRecord const & record) const {
754  return record.get(_psfShapeResult);
755 }
756 
758  record.set(_shapeResult, value);
759  record.set(_centroidResult, value);
760  record.set(_fluxResult, value);
761  record.set(_flux_xx_Cov, value.flux_xx_Cov);
762  record.set(_flux_yy_Cov, value.flux_yy_Cov);
763  record.set(_flux_xy_Cov, value.flux_xy_Cov);
764  for (size_t n = 0; n < SdssShapeAlgorithm::N_FLAGS; ++n) {
765  if (n == SdssShapeAlgorithm::PSF_SHAPE_BAD.number && !_includePsf) continue;
766  _flagHandler.setValue(record, n, value.flags[n]);
767  }
768 }
769 
771  afwGeom::ellipses::Quadrupole const & value) const {
772  record.set(_psfShapeResult, value);
773 }
774 
776  return _shapeResult == other._shapeResult &&
777  _centroidResult == other._centroidResult &&
778  _fluxResult == other._fluxResult &&
779  _psfShapeResult == other._psfShapeResult &&
780  _flux_xx_Cov == other._flux_xx_Cov &&
781  _flux_yy_Cov == other._flux_yy_Cov &&
782  _flux_xy_Cov == other._flux_xy_Cov;
783  // don't bother with flags - if we've gotten this far, it's basically impossible the flags don't match
784 }
785 
787  return _shapeResult.isValid() &&
788  _centroidResult.isValid() &&
789  _fluxResult.isValid() &&
790  _psfShapeResult.isValid() &&
791  _flux_xx_Cov.isValid() &&
792  _flux_yy_Cov.isValid() &&
793  _flux_xy_Cov.isValid();
794  // don't bother with flags - if we've gotten this far, it's basically impossible the flags are invalid
795 }
796 
798  Control const & ctrl,
799  std::string const & name,
800  afw::table::Schema & schema
801 )
802  : _ctrl(ctrl),
803  _resultKey(ResultKey::addFields(schema, name, ctrl.doMeasurePsf)),
804  _centroidExtractor(schema, name)
805 {}
806 
807 template <typename ImageT>
809  ImageT const & image,
810  afw::geom::Point2D const & center,
811  bool negative,
812  Control const & control
813 ) {
814  double xcen = center.getX(); // object's column position
815  double ycen = center.getY(); // object's row position
816 
817  xcen -= image.getX0(); // work in image Pixel coordinates
818  ycen -= image.getY0();
819 
820  float shiftmax = control.maxShift; // Max allowed centroid shift
821  if (shiftmax < 2) {
822  shiftmax = 2;
823  } else if (shiftmax > 10) {
824  shiftmax = 10;
825  }
826 
827  SdssShapeResult result;
828  try {
829  result.flags[FAILURE.number] = !getAdaptiveMoments(
830  image, control.background, xcen, ycen, shiftmax, &result,
831  control.maxIter, control.tol1, control.tol2, negative
832  );
833  } catch (pex::exceptions::Exception & err) {
834  result.flags[FAILURE.number] = true;
835  }
836  if (result.flags[UNWEIGHTED.number] || result.flags[SHIFT.number]) {
837  // These are also considered fatal errors in terms of the quality of the results,
838  // even though they do produce some results.
839  result.flags[FAILURE.number] = true;
840  }
841  if (result.getQuadrupole().getIxx()*result.getQuadrupole().getIyy() <
842  (1.0 + 1.0e-6)*result.getQuadrupole().getIxy()*result.getQuadrupole().getIxy())
843  // We are checking that Ixx*Iyy > (1 + epsilon)*Ixy*Ixy where epsilon is suitably small. The
844  // value of epsilon used here is a magic number. DM-5801 is supposed to figure out if we are
845  // to keep this value.
846  {
847  if (!result.flags[FAILURE.number]) {
848  throw LSST_EXCEPT(
850  "Should not get singular moments unless a flag is set");
851  }
852  }
853 
854  // getAdaptiveMoments() just computes the zeroth moment in result.flux (and its error in
855  // result.fluxSigma, result.flux_xx_Cov, etc.) That's related to the flux by some geometric
856  // factors, which we apply here.
857  double fluxScale = computeFluxScale(result);
858 
859  result.flux *= fluxScale;
860  result.fluxSigma *= fluxScale;
861  result.x += image.getX0();
862  result.y += image.getY0();
863 
864  if (ImageAdaptor<ImageT>::hasVariance) {
865  result.flux_xx_Cov *= fluxScale;
866  result.flux_yy_Cov *= fluxScale;
867  result.flux_xy_Cov *= fluxScale;
868  }
869 
870  return result;
871 }
872 
873 template <typename ImageT>
875  ImageT const & image,
876  afw::geom::ellipses::Quadrupole const & shape,
877  afw::geom::Point2D const & center
878 ) {
879  // while arguments to computeFixedMomentsFlux are in PARENT coordinates, the implementation is LOCAL.
880  afw::geom::Point2D localCenter = center - afw::geom::Extent2D(image.getXY0());
881 
882  afwGeom::BoxI const bbox = computeAdaptiveMomentsBBox(image.getBBox(afw::image::LOCAL),
883  localCenter,
884  shape.getIxx(), shape.getIxy(), shape.getIyy());
885 
886  std::tuple<std::pair<bool, double>, double, double, double> weights =
887  getWeights(shape.getIxx(), shape.getIxy(), shape.getIyy());
888 
889  FluxResult result;
890 
891  if (!std::get<0>(weights).first) {
892  throw pex::exceptions::InvalidParameterError("Input shape is singular");
893  }
894 
895  double const w11 = std::get<1>(weights);
896  double const w12 = std::get<2>(weights);
897  double const w22 = std::get<3>(weights);
898  bool const interp = shouldInterp(shape.getIxx(), shape.getIyy(), std::get<0>(weights).second);
899 
900  double i0 = 0; // amplitude of Gaussian
901  if (calcmom<true>(ImageAdaptor<ImageT>().getImage(image), localCenter.getX(), localCenter.getY(),
902  bbox, 0.0, interp, w11, w12, w22, &i0, NULL, NULL, NULL, NULL, NULL, NULL, NULL)< 0) {
903  throw LSST_EXCEPT(pex::exceptions::RuntimeError, "Error from calcmom");
904  }
905 
906  double const wArea = afw::geom::PI*std::sqrt(shape.getDeterminant());
907 
908  result.flux = i0*2*wArea;
909 
910  if (ImageAdaptor<ImageT>::hasVariance) {
911  int ix = static_cast<int>(center.getX() - image.getX0());
912  int iy = static_cast<int>(center.getY() - image.getY0());
913  if (!image.getBBox(afw::image::LOCAL).contains(afw::geom::Point2I(ix, iy))) {
915  (boost::format("Center (%d,%d) not in image (%dx%d)") %
916  ix % iy % image.getWidth() % image.getHeight()).str());
917  }
918  double var = ImageAdaptor<ImageT>().getVariance(image, ix, iy);
919  double i0Err = std::sqrt(var/wArea);
920  result.fluxSigma = i0Err*2*wArea;
921  }
922 
923  return result;
924 }
925 
927  afw::table::SourceRecord & measRecord,
928  afw::image::Exposure<float> const & exposure
929 ) const {
930  bool negative = false;
931 
932  try {
933  negative = measRecord.get(measRecord.getSchema().find<afw::table::Flag>("flags_negative").key);
934  } catch(pexExcept::Exception &e) {
935  }
937  exposure.getMaskedImage(),
938  _centroidExtractor(measRecord, _resultKey.getFlagHandler()),
939  negative,
940  _ctrl
941  );
942 
943  if (_ctrl.doMeasurePsf) {
944  // Compute moments of Psf model. In the interest of implementing this quickly, we're just
945  // calling Psf::computeShape(), which delegates to SdssShapeResult::computeAdaptiveMoments
946  // for all nontrivial Psf classes. But this could in theory save the results of a shape
947  // computed some other way as part of base_SdssShape, which might be confusing. We should
948  // fix this eventually either by making Psf shape measurement not part of base_SdssShape, or
949  // by making the measurements stored with shape.sdss always computed via the
950  // SdssShapeAlgorithm instead of delegating to the Psf class.
951  try {
952  PTR(afw::detection::Psf const) psf = exposure.getPsf();
953  if (!psf) {
954  result.flags[PSF_SHAPE_BAD.number] = true;
955  } else {
956  _resultKey.setPsfShape(measRecord, psf->computeShape(afw::geom::Point2D(result.x, result.y)));
957  }
958  } catch (pex::exceptions::Exception & err) {
959  result.flags[PSF_SHAPE_BAD.number] = true;
960  }
961  }
962 
963  measRecord.set(_resultKey, result);
964 }
965 
967  afw::table::SourceRecord & measRecord,
968  MeasurementError * error
969 ) const {
970  _resultKey.getFlagHandler().handleFailure(measRecord, error);
971 }
972 
973 #define INSTANTIATE_IMAGE(IMAGE) \
974  template SdssShapeResult SdssShapeAlgorithm::computeAdaptiveMoments( \
975  IMAGE const &, \
976  afw::geom::Point2D const &, \
977  bool, \
978  Control const & \
979  ); \
980  template FluxResult SdssShapeAlgorithm::computeFixedMomentsFlux( \
981  IMAGE const &, \
982  afw::geom::ellipses::Quadrupole const &, \
983  afw::geom::Point2D const & \
984  )
985 
986 #define INSTANTIATE_PIXEL(PIXEL) \
987  INSTANTIATE_IMAGE(lsst::afw::image::Image<PIXEL>); \
988  INSTANTIATE_IMAGE(lsst::afw::image::MaskedImage<PIXEL>);
989 
990 INSTANTIATE_PIXEL(int);
991 INSTANTIATE_PIXEL(float);
992 INSTANTIATE_PIXEL(double);
993 
995  Control const & ctrl,
996  std::string const & name,
997  afw::table::SchemaMapper & mapper
998 ) :
999  BaseTransform{name},
1000  _fluxTransform{name, mapper},
1001  _centroidTransform{name, mapper}
1002 {
1003  // If the input schema has a PSF flag -- it's optional -- assume we are also transforming the PSF.
1004  _transformPsf = mapper.getInputSchema().getNames().count("sdssShape_flag_psf") ? true : false;
1005 
1006  // Skip the last flag if not transforming the PSF shape.
1007  for (std::size_t i = 0; i < SdssShapeAlgorithm::getFlagDefinitions().size(); i++) {
1009  if (flag == SdssShapeAlgorithm::FAILURE) continue;
1010  if (mapper.getInputSchema().getNames().count(
1011  mapper.getInputSchema().join(name, flag.name)) == 0) continue;
1012  afw::table::Key<afw::table::Flag> key = mapper.getInputSchema().find<afw::table::Flag>(
1013  name + "_" + flag.name).key;
1014  mapper.addMapping(key);
1015  }
1016 
1017  _outShapeKey = ShapeResultKey::addFields(mapper.editOutputSchema(), name, "Shape in celestial moments",
1019  if (_transformPsf) {
1020  _outPsfShapeKey = afwTable::QuadrupoleKey::addFields(mapper.editOutputSchema(), name + "_psf",
1021  "PSF shape in celestial moments",
1023  }
1024 }
1025 
1027  afw::table::SourceCatalog const & inputCatalog,
1028  afw::table::BaseCatalog & outputCatalog,
1029  afw::image::Wcs const & wcs,
1030  afw::image::Calib const & calib
1031 ) const {
1032  // The flux and cetroid transforms will check that the catalog lengths
1033  // match and throw if not, so we don't repeat the test here.
1034  _fluxTransform(inputCatalog, outputCatalog, wcs, calib);
1035  _centroidTransform(inputCatalog, outputCatalog, wcs, calib);
1036 
1037  CentroidResultKey centroidKey(inputCatalog.getSchema()[_name]);
1038  ShapeResultKey inShapeKey(inputCatalog.getSchema()[_name]);
1039  afwTable::QuadrupoleKey inPsfShapeKey;
1040  if (_transformPsf) {
1041  inPsfShapeKey = afwTable::QuadrupoleKey(
1042  inputCatalog.getSchema()[inputCatalog.getSchema().join(_name, "psf")]);
1043  }
1044 
1045  afw::table::SourceCatalog::const_iterator inSrc = inputCatalog.begin();
1046  afw::table::BaseCatalog::iterator outSrc = outputCatalog.begin();
1047  for (; inSrc != inputCatalog.end(); ++inSrc, ++outSrc) {
1048  ShapeResult inShape = inShapeKey.get(*inSrc);
1049  ShapeResult outShape;
1050 
1051  // The transformation from the (x, y) to the (Ra, Dec) basis.
1052  afw::geom::AffineTransform crdTr = wcs.linearizePixelToSky(centroidKey.get(*inSrc).getCentroid(),
1054  outShape.setShape(inShape.getShape().transform(crdTr.getLinear()));
1055 
1056  // Transformation matrix from pixel to celestial basis.
1058  outShape.setShapeErr((m*inShape.getShapeErr().cast<double>()*m.transpose()).cast<ErrElement>());
1059 
1060  _outShapeKey.set(*outSrc, outShape);
1061 
1062  if (_transformPsf) {
1063  _outPsfShapeKey.set(*outSrc, inPsfShapeKey.get(*inSrc).transform(crdTr.getLinear()));
1064  }
1065  }
1066 }
1067 
1068 }}} // end namespace lsst::meas::base
int y
int iter
std::size_t size() const
return the current size (number of defined elements) of the collection
Definition: FlagHandler.h:133
virtual void set(afw::table::BaseRecord &record, ShapeResult const &value) const
Set a ShapeResult in the given record.
bool doMeasurePsf
"Whether to also compute the shape of the PSF model" ;
Definition: SdssShape.h:57
virtual void measure(afw::table::SourceRecord &measRecord, afw::image::Exposure< float > const &exposure) const
Called to measure a single child source in an image.
Definition: SdssShape.cc:926
static Result computeAdaptiveMoments(ImageT const &image, afw::geom::Point2D const &position, bool negative=false, Control const &ctrl=Control())
Compute the adaptive Gaussian-weighted moments of an image.
ShapeTrMatrix makeShapeTransformMatrix(afw::geom::LinearTransform const &xform)
Construct a matrix suitable for transforming second moments.
virtual void set(BaseRecord &record, geom::ellipses::Quadrupole const &value) const
afw::table::Key< afw::table::Array< ImagePixelT > > image
ErrElement yy_xy_Cov
yy,xy term in the uncertainty convariance matrix
Extent< double, 2 > Extent2D
static FluxResult computeFixedMomentsFlux(ImageT const &image, afw::geom::ellipses::Quadrupole const &shape, afw::geom::Point2D const &position)
Compute the flux within a fixed Gaussian aperture.
Definition: SdssShape.cc:874
float tol2
"Convergence tolerance for FWHM" ;
Definition: SdssShape.h:56
int positionToIndex(double pos)
static FlagDefinition const SHIFT
Definition: SdssShape.h:161
T exp(T... args)
AngleUnit constexpr radians
std::string join(std::string const &a, std::string const &b) const
Shape const getShape() const
Return an afw::geom::ellipses object corresponding to xx, yy, xy.
iterator at(int const x, int const y) const
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
LinearTransform const & getLinear() const
A reusable struct for centroid measurements.
X
afw::geom::ellipses::Quadrupole getQuadrupole()
double background
"Additional value to add to background" ;
Definition: SdssShape.h:52
bool isValid() const
Return True if the centroid key is valid.
geom::AffineTransform linearizePixelToSky(coord::Coord const &coord, geom::AngleUnit skyUnit=geom::degrees) const
The full covariance matrix is provided.
Definition: constants.h:44
STL namespace.
CentroidElement x
x (column) coordinate of the measured position
Y
T make_tuple(T... args)
#define PTR(...)
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...
Definition: SdssShape.cc:966
int maxIter
"Maximum number of iterations" ;
Definition: SdssShape.h:53
bool operator==(SdssShapeResultKey const &other) const
Compare the FunctorKey for equality with another, using the underlying Keys.
Definition: SdssShape.cc:775
SchemaItem< T > find(std::string const &name) const
Transformer transform(LinearTransform const &transform)
void setValue(afw::table::BaseRecord &record, std::size_t i, bool value) const
Set the flag field corresponding to the given flag index.
Definition: FlagHandler.h:276
Key< T > addField(Field< T > const &field, bool doReplace=false)
Exception to be thrown when a measurement algorithm experiences a known failure mode.
Definition: exceptions.h:48
static FlagDefinition const UNWEIGHTED
Definition: SdssShape.h:160
ErrElement xySigma
1-Sigma uncertainty on xy (sqrt of variance)
static FluxResultKey addFields(afw::table::Schema &schema, std::string const &name, std::string const &doc)
Add a pair of _flux, _fluxSigma fields to a Schema, and return a FluxResultKey that points to them...
Field< T >::Value get(Key< T > const &key) const
float tol1
"Convergence tolerance for e1,e2" ;
Definition: SdssShape.h:55
bool isValid() const
Return True if the key is valid.
Definition: SdssShape.cc:786
static FlagDefinition const UNWEIGHTED_BAD
Definition: SdssShape.h:159
string name
double maxShift
"Maximum centroid shift, limited to 2-10" ;
Definition: SdssShape.h:54
float ErrElement
Definition: constants.h:53
ErrElement flux_yy_Cov
flux, yy term in the uncertainty covariance matrix
Definition: SdssShape.h:241
bool isValid() const
A FunctorKey for ShapeResult.
STL class.
T min(T... args)
double constexpr PI
ErrElement xx_xy_Cov
xx,xy term in the uncertainty convariance matrix
void setShapeErr(ShapeCov const &matrix)
Set the struct uncertainty elements from the given matrix, with rows and columns ordered (xx...
ErrElement xx_yy_Cov
xx,yy term in the uncertainty convariance matrix
Utility class for handling flag fields that indicate the failure modes of an algorithm.
Definition: FlagHandler.h:156
A reusable struct for moments-based shape measurements.
virtual void operator()(afw::table::SourceCatalog const &inputCatalog, afw::table::BaseCatalog &outputCatalog, afw::image::Wcs const &wcs, afw::image::Calib const &calib) const
Definition: SdssShape.cc:1026
first
static unsigned int const N_FLAGS
Definition: SdssShape.h:157
MaskedImageT getMaskedImage()
static QuadrupoleKey addFields(Schema &schema, std::string const &name, std::string const &doc, CoordinateType coordType=CoordinateType::PIXEL)
Result object SdssShapeAlgorithm.
Definition: SdssShape.h:238
Point< double, 2 > Point2D
static FlagDefinitionList const & getFlagDefinitions()
Definition: SdssShape.cc:57
SdssShapeResult()
Constructor; initializes everything to NaN.
Definition: SdssShape.cc:660
SdssShapeAlgorithm(Control const &ctrl, std::string const &name, afw::table::Schema &schema)
Definition: SdssShape.cc:797
T make_pair(T... args)
Flux flux
Measured flux in DN.
Definition: FluxUtilities.h:38
virtual void set(afw::table::BaseRecord &record, SdssShapeResult const &value) const
Set an SdssShapeResult in the given record.
Definition: SdssShape.cc:757
static CentroidResultKey addFields(afw::table::Schema &schema, std::string const &name, std::string const &doc, UncertaintyEnum uncertainty)
Add the appropriate fields to a Schema, and return a CentroidResultKey that manages them...
A C++ control class to handle SdssShapeAlgorithm&#39;s configuration.
Definition: SdssShape.h:50
T max(T... args)
bool isValid() const
Return True if the shape key is valid.
bool isValid() const
Return True if both the flux and fluxSigma Keys are valid.
static FlagHandler addFields(afw::table::Schema &schema, std::string const &prefix, FlagDefinitionList const &flagDefs, FlagDefinitionList const &exclDefs=FlagDefinitionList::getEmptyList())
Add Flag fields to a schema, creating a FlagHandler object to manage them.
Definition: FlagHandler.cc:37
int x
std::bitset< SdssShapeAlgorithm::N_FLAGS > flags
Status flags (see SdssShapeAlgorithm).
Definition: SdssShape.h:244
virtual SdssShapeResult get(afw::table::BaseRecord const &record) const
Get an SdssShapeResult from the given record.
Definition: SdssShape.cc:738
virtual afw::geom::ellipses::Quadrupole getPsfShape(afw::table::BaseRecord const &record) const
Get a Quadrupole for the Psf from the given record.
Definition: SdssShape.cc:753
second
#define LSST_EXCEPT(type,...)
Base::const_iterator const_iterator
static FlagDefinition const MAXITER
Definition: SdssShape.h:162
std::shared_ptr< lsst::afw::detection::Psf > getPsf()
T pow(T... args)
static FlagDefinition const PSF_SHAPE_BAD
Definition: SdssShape.h:163
A FunctorKey that maps SdssShapeResult to afw::table Records.
Definition: SdssShape.h:71
Algorithm provides no uncertainy information at all.
Definition: constants.h:42
Abstract base class for all C++ measurement transformations.
Definition: Transform.h:82
static FlagDefinition const FAILURE
Definition: SdssShape.h:158
void handleFailure(afw::table::BaseRecord &record, MeasurementError const *error=nullptr) const
Handle an expected or unexpected Exception thrown by a measurement algorithm.
Definition: FlagHandler.cc:93
ErrElement yySigma
1-Sigma uncertainty on yy (sqrt of variance)
A FunctorKey for CentroidResult.
T isnan(T... args)
void set(Key< T > const &key, U const &value)
table::Box2IKey bbox
CentroidElement y
y (row) coordinate of the measured position
ErrElement xxSigma
1-Sigma uncertainty on xx (sqrt of variance)
void clip(Box2I const &other)
bool getValue(afw::table::BaseRecord const &record, std::size_t i) const
Return the value of the flag field corresponding to the given flag index.
Definition: FlagHandler.h:256
T quiet_NaN(T... args)
T sqrt(T... args)
m
SdssShapeResultKey()
Default constructor; instance will not be usuable unless subsequently assigned to.
Definition: SdssShape.h:92
FlagHandler const & getFlagHandler() const
Definition: SdssShape.h:127
#define INSTANTIATE_PIXEL(PIXEL)
Definition: SdssShape.cc:986
vector-type utility class to build a collection of FlagDefinitions
Definition: FlagHandler.h:63
ErrElement flux_xx_Cov
flux, xx term in the uncertainty covariance matrix
Definition: SdssShape.h:240
void setShape(Shape const &shape)
Set struct elements from the given Quadrupole object.
FluxErrElement fluxSigma
1-Sigma error (sqrt of variance) on flux in DN.
Definition: FluxUtilities.h:39
Eigen::Matrix< ShapeElement, 3, 3, Eigen::DontAlign > ShapeTrMatrix
Definition: constants.h:60
static SdssShapeResultKey addFields(afw::table::Schema &schema, std::string const &name, bool doMeasurePsf)
Add the appropriate fields to a Schema, and return a SdssShapeResultKey that manages them...
Definition: SdssShape.cc:667
virtual void setPsfShape(afw::table::BaseRecord &record, afw::geom::ellipses::Quadrupole const &value) const
Set a Quadrupole for the Psf at the position of the given record.
Definition: SdssShape.cc:770
SdssShapeTransform(Control const &ctrl, std::string const &name, afw::table::SchemaMapper &mapper)
Definition: SdssShape.cc:994
A reusable result struct for flux measurements.
Definition: FluxUtilities.h:37
double constexpr TWOPI
ErrElement flux_xy_Cov
flux, xy term in the uncertainty covariance matrix
Definition: SdssShape.h:242
ShapeCov const getShapeErr() const
Return the 3x3 symmetric covariance matrix, with rows and columns ordered (xx, yy, xy)
static ShapeResultKey addFields(afw::table::Schema &schema, std::string const &name, std::string const &doc, UncertaintyEnum uncertainty, afw::table::CoordinateType coordType=afw::table::CoordinateType::PIXEL)
Add the appropriate fields to a Schema, and return a ShapeResultKey that manages them.