26#include "boost/math/special_functions/bessel.hpp"
27#include "boost/shared_array.hpp"
42inline double J1(
double const x) {
return boost::math::cyl_bessel_j(1,
x); }
46inline T sinc(T
const x) {
53template <
typename CoordT>
54class CircularAperture {
59 CoordT
const taperwidth
63 _taperwidth1(taperwidth),
64 _taperwidth2(taperwidth),
65 _k1(1.0 / (2.0 * taperwidth)),
66 _k2(1.0 / (2.0 * taperwidth)),
67 _taperLo1(radius1 - 0.5 * taperwidth),
68 _taperHi1(radius1 + 0.5 * taperwidth),
69 _taperLo2(radius2 - 0.5 * taperwidth),
70 _taperHi2(radius2 + 0.5 * taperwidth) {
75 if (_radius1 > _radius2) {
77 pex::exceptions::InvalidParameterError,
78 (boost::format(
"rad2 less than rad1: (rad1=%.2f, rad2=%.2f) ") % _radius1 % _radius2)
81 if (_radius1 < 0.0 || _radius2 < 0.0) {
83 pex::exceptions::InvalidParameterError,
84 (boost::format(
"radii must be >= 0 (rad1=%.2f, rad2=%.2f) ") % _radius1 % _radius2)
94 if (_radius1 < 0.5 * _taperwidth1) {
95 _taperwidth1 = 2.0 * _radius1;
96 _k1 = 1.0 / (2.0 * _taperwidth1);
100 if ((_radius2 - _radius1) < 0.5 * (_taperwidth1 + _taperwidth2)) {
103 if ((_radius2 - _radius2) < 0.5 * _taperwidth1) {
104 _taperwidth1 = _taperwidth2 = 0.5 * (_radius2 - _radius1);
105 _k1 = _k2 = 1.0 / (2.0 * _taperwidth1);
109 _taperwidth2 = _radius2 - _radius1 - _taperwidth1;
110 _k2 = 1.0 / (2.0 * _taperwidth2);
113 _taperLo1 = _radius1 - 0.5 * _taperwidth1;
114 _taperHi1 = _radius1 + 0.5 * _taperwidth1;
115 _taperLo2 = _radius2 - 0.5 * _taperwidth2;
116 _taperHi2 = _radius2 + 0.5 * _taperwidth2;
122 CoordT operator()(CoordT
const x, CoordT
const y)
const {
124 if (xyrad < _taperLo1) {
126 }
else if (xyrad >= _taperLo1 && xyrad <= _taperHi1) {
128 }
else if (xyrad > _taperHi1 && xyrad <= _taperLo2) {
130 }
else if (xyrad > _taperLo2 && xyrad <= _taperHi2) {
137 CoordT getRadius1() {
return _radius1; }
138 CoordT getRadius2() {
return _radius2; }
141 CoordT _radius1, _radius2;
142 CoordT _taperwidth1, _taperwidth2;
144 CoordT _taperLo1, _taperHi1;
145 CoordT _taperLo2, _taperHi2;
148template <
typename CoordT>
151 CircApPolar(
double radius,
double taperwidth) : _ap(CircularAperture<CoordT>(0.0,
radius, taperwidth)) {}
152 CoordT operator()(
double r)
const {
return r * _ap(r, 0.0); }
155 CircularAperture<CoordT> _ap;
161template <
typename IntegrandT>
164 SincAperture(CircularAperture<IntegrandT>
const& ap,
168 : _ap(ap), _ix(ix), _iy(iy) {}
170 IntegrandT operator()(IntegrandT
const x, IntegrandT
const y)
const {
171 double const fourierConvention =
geom::PI;
172 double const dx = fourierConvention * (
x - _ix);
173 double const dy = fourierConvention * (
y - _iy);
174 double const fx = sinc<double>(dx);
175 double const fy = sinc<double>(dy);
176 return (1.0 + _ap(
x,
y) * fx * fy);
180 CircularAperture<IntegrandT>
const& _ap;
184class GaussPowerFunctor {
186 GaussPowerFunctor(
double sigma) : _sigma(
sigma) {}
188 double operator()(
double kx,
double ky)
const {
189 double const k = ::sqrt(kx * kx + ky * ky);
190 double const gauss = ::exp(-0.5 * k * k * _sigma * _sigma);
191 return gauss * gauss;
199 GaussPowerFunctor gaussPower(
sigma);
208 double trueLeak = (powerInf - truePower) / powerInf;
218template <
typename PixelT>
220 double const taperwidth) {
223 int log2 =
static_cast<int>(::ceil(::log10(2.0 * rad2) /
log10(2.0)));
227 int hwid =
pow(2, log2);
228 int width = 2 * hwid - 1;
230 int const xwidth = width;
231 int const ywidth = width;
233 int const x0 = -xwidth / 2;
234 int const y0 = -ywidth / 2;
237 auto coeffImage = std::make_shared<afw::image::Image<PixelT>>(
geom::ExtentI(xwidth, ywidth), initweight);
238 coeffImage->setXY0(x0, y0);
242 double tolerance = 1.0e-12;
244 double err = 2.0 * tolerance;
245 double apEff =
geom::PI * rad2 * rad2;
249 while (err > tolerance && i < maxIt) {
250 CircApPolar<double> apPolar1(radIn, taperwidth);
251 CircApPolar<double> apPolar2(radIn + dr, taperwidth);
254 double dadr = (a2 - a1) / dr;
255 double radNew = radIn - (a1 - apEff) / dadr;
256 err = (a1 - apEff) / apEff;
260 CircularAperture<double> ap(rad1, rad2, taperwidth);
266 double const limit = rad2 + taperwidth;
267 double const x1 = -limit;
268 double const x2 = limit;
269 double const y1 = -limit;
270 double const y2 = limit;
272 for (
int iY = y0; iY != y0 + coeffImage->getHeight(); ++iY) {
278 SincAperture<double> sincAp(ap, iX, iY);
285 double const dx = iX;
286 double const dy = iY;
287 *
ptr = (
std::sqrt(dx * dx + dy * dy) < xwidth / 2) ? integral - (x2 - x1) * (y2 - y1) : 0.0;
293 for (
int iY = y0; iY != y0 + coeffImage->getHeight(); ++iY) {
302 coeffImage->writeFits(
"cimage.fits");
310 FftShifter(
int xwid) : _xwid(xwid) {}
312 if (
x >= _xwid / 2) {
313 return x - _xwid / 2;
315 return x + _xwid / 2 + 1;
324 double c = ::cos(
angle);
325 double s = ::sin(
angle);
335template <
typename PixelT>
338 double const ellipticity) {
341 int log2 =
static_cast<int>(::ceil(::log10(2.0 * rad2) /
log10(2.0)));
345 int hwid =
pow(2, log2);
346 int wid = 2 * hwid - 1;
347 int xcen = wid / 2, ycen = wid / 2;
348 FftShifter fftshift(wid);
354 fftw_plan plan = fftw_plan_dft_2d(wid, wid,
reinterpret_cast<fftw_complex*
>(c),
355 reinterpret_cast<fftw_complex*
>(c), FFTW_BACKWARD, FFTW_ESTIMATE);
360 double const scale = (1.0 - ellipticity);
361 for (
int iY = 0; iY < wid; ++iY) {
362 int const fY = fftshift.shift(iY);
363 double const ky = (
static_cast<double>(iY) - ycen) / wid;
365 for (
int iX = 0; iX < wid; ++iX) {
366 int const fX = fftshift.shift(iX);
367 double const kx =
static_cast<double>(iX - xcen) / wid;
371 double kxr = coo.first;
372 double kyr = coo.second;
374 double const k = ::sqrt(kxr * kxr + scale * scale * kyr * kyr);
376 double const airy1 = (rad1 > 0 ? rad1 * J1(twoPiRad1 * k) : 0.0) / k;
377 double const airy2 = rad2 * J1(twoPiRad2 * k) / k;
378 double const airy = airy2 - airy1;
387 fftw_destroy_plan(plan);
390 auto coeffImage = std::make_shared<afw::image::Image<PixelT>>(
geom::ExtentI(wid, wid), 0.0);
392 for (
int iY = 0; iY != coeffImage->getHeight(); ++iY) {
397 int fX = fftshift.shift(iX);
398 int fY = fftshift.shift(iY);
399 *
ptr =
static_cast<PixelT>(c[fY * wid + fX].
real() / (wid * wid));
405 coeffImage->setXY0(-wid / 2, -wid / 2);
413template <
typename PixelT>
417 int log2 =
static_cast<int>(::ceil(::log10(2.0 * rad2) /
log10(2.0)));
421 int hwid =
pow(2, log2);
422 int wid = 2 * hwid - 1;
423 int xcen = wid / 2, ycen = wid / 2;
424 FftShifter fftshift(wid);
431 fftw_plan plan = fftw_plan_dft_2d(wid, wid,
reinterpret_cast<fftw_complex*
>(c),
432 reinterpret_cast<fftw_complex*
>(c), FFTW_BACKWARD, FFTW_ESTIMATE);
436 for (
int iY = 0; iY < wid; ++iY) {
437 int const fY = fftshift.shift(iY);
438 double const ky = (
static_cast<double>(iY) - ycen) / wid;
440 for (
int iX = 0; iX < wid; ++iX) {
441 int const fX = fftshift.shift(iX);
444 double const iXcen =
static_cast<double>(iX - xcen);
445 double const kx = iXcen / wid;
447 double const k = ::sqrt(kx * kx + ky * ky);
448 double const airy1 = (rad1 > 0 ? rad1 * J1(twoPiRad1 * k) : 0.0) / k;
449 double const airy2 = rad2 * J1(twoPiRad2 * k) / k;
450 double const airy = airy2 - airy1;
451 c[fY * wid + fX] = airy;
454 int fxy = fftshift.shift(wid / 2);
455 c[fxy * wid + fxy] = {
geom::PI * (rad2 * rad2 - rad1 * rad1), 0.};
459 fftw_destroy_plan(plan);
462 auto coeffImage = std::make_shared<afw::image::Image<PixelT>>(
geom::ExtentI(wid, wid), 0.0);
464 for (
int iY = 0; iY != coeffImage->getHeight(); ++iY) {
470 int fX = iX < hwid ? hwid - iX - 1 : iX - hwid + 1;
471 int fY = iY < hwid ? hwid - iY - 1 : iY - hwid + 1;
472 *
ptr =
static_cast<PixelT>(c[fY * wid + fX].
real() / (wid * wid));
478 coeffImage->setXY0(-wid / 2, -wid / 2);
484template <
typename PixelT>
485SincCoeffs<PixelT>& SincCoeffs<PixelT>::getInstance() {
486 static SincCoeffs<PixelT> instance;
490template <
typename PixelT>
492 if (r1 < 0.0 || r2 < r1) {
494 (boost::format(
"Invalid r1,r2 = %f,%f") % r1 % r2).
str());
496 double const innerFactor = r1 / r2;
498 if (!getInstance()._lookup(axes, innerFactor)) {
500 getInstance()._cache[r2][innerFactor] =
coeff;
504template <
typename PixelT>
508 return coeff ?
coeff : calculate(axes, innerFactor);
511template <
typename PixelT>
514 if (innerFactor < 0.0 || innerFactor > 1.0) {
516 (boost::format(
"innerFactor = %f is not between 0 and 1") % innerFactor).
str());
522 if (!FuzzyCompare<float>().isEqual(axes.
getA(), axes.
getB())) {
525 typename CoeffMapMap::const_iterator iter1 = _cache.find(axes.
getA());
526 if (iter1 == _cache.end()) {
529 typename CoeffMap::const_iterator iter2 = iter1->second.find(innerFactor);
530 return (iter2 == iter1->second.end()) ? null : iter2->second;
533template <
typename PixelT>
536 if (innerFactor < 0.0 || innerFactor > 1.0) {
538 (boost::format(
"innerFactor = %f is not between 0 and 1") % innerFactor).
str());
544 double const rad1 = axes.
getA() * innerFactor;
545 double const rad2 = axes.
getA();
547 if (FuzzyCompare<float>().isEqual(axes.
getA(), axes.
getB())) {
549 return calcImageKSpaceReal<PixelT>(rad1, rad2);
552 double const ellipticity = 1.0 - axes.
getB() / axes.
getA();
553 return calcImageKSpaceCplx<PixelT>(rad1, rad2, axes.
getTheta(), ellipticity);
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
table::Key< double > angle
afw::table::Key< double > sigma
An ellipse core for the semimajor/semiminor axis and position angle parametrization (a,...
double const getTheta() const
double const getA() const
double const getB() const
typename _view_t::x_iterator x_iterator
An iterator for traversing the pixels in a row.
A singleton to calculate and cache the coefficients for sinc photometry.
static std::shared_ptr< CoeffT > calculate(afw::geom::ellipses::Axes const &outerEllipse, double const innerFactor=0.0)
Calculate the coefficients for an aperture.
static void cache(float rInner, float rOuter)
Cache the coefficients for a particular aperture.
static std::shared_ptr< CoeffT const > get(afw::geom::ellipses::Axes const &outerEllipse, float const innerRadiusFactor=0.0)
Get the coefficients for an aperture.
Reports invalid arguments.
def scale(algorithm, min, max=None, frame=None)
auto integrate(UnaryFunctionT func, Arg const a, Arg const b, double eps=1.0e-6)
The 1D integrator.
auto integrate2d(BinaryFunctionT func, X x1, X x2, Y y1, Y y2, double eps=1.0e-6)
The 2D integrator.
double constexpr PI
The ratio of a circle's circumference to diameter.
table::Key< table::Array< double > > coeff