lsst.meas.extensions.psfex  15.0-6-g86a76e2+8
PsfexPsf.cc
Go to the documentation of this file.
1 // -*- LSST-C++ -*-
2 
3 /*
4  * LSST Data Management System
5  * Copyright 2008, 2009, 2010 LSST Corporation.
6  *
7  * This product includes software developed by the
8  * LSST Project (http://www.lsst.org/).
9  *
10  * This program is free software: you can redistribute it and/or modify
11  * it under the terms of the GNU General Public License as published by
12  * the Free Software Foundation, either version 3 of the License, or
13  * (at your option) any later version.
14  *
15  * This program is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the LSST License Statement and
21  * the GNU General Public License along with this program. If not,
22  * see <http://www.lsstcorp.org/LegalNotices/>.
23  */
24 
30 #include <cmath>
31 #include <cassert>
32 #include <numeric>
33 #include <memory>
34 
35 #include "lsst/base.h"
36 #include "lsst/pex/exceptions.h"
39 extern "C" {
40 #include "vignet.h"
41 }
45 
46 namespace lsst { namespace meas { namespace extensions { namespace psfex {
47 
48 namespace afw = lsst::afw;
49 
50 PsfexPsf::PsfexPsf(
51  lsst::meas::extensions::psfex::Psf const& psf,
52  afw::geom::Point2D const & averagePosition
53  ) : ImagePsf(), _averagePosition(averagePosition),
54  _size(psf.impl->dim),
55  _comp(psf.impl->npix),
56  _context(psf.impl->poly->ndim)
57 
58 {
59  _poly = poly_copy(psf.impl->poly);
60 
61  _pixstep = psf.impl->pixstep;
62 
63  std::copy(psf.impl->size, psf.impl->size + psf.impl->dim, _size.begin());
64 
65  std::copy(psf.impl->comp, psf.impl->comp + psf.impl->npix, _comp.begin());
66 
67  for (int i = 0; i != psf.impl->poly->ndim; ++i) {
68  _context[i].first = psf.impl->contextoffset[i];
69  _context[i].second = psf.impl->contextscale[i];
70  }
71 }
72 
73  PsfexPsf::PsfexPsf() : ImagePsf(),
74  _averagePosition(afw::geom::Point2I(0, 0)),
75  _poly(0),
76  _pixstep(0.0),
77  _size(),
78  _comp(),
79  _context()
80 {
81  ;
82 }
83 
85 {
86  poly_end(_poly);
87 }
88 
90 PsfexPsf::clone() const {
91  return std::make_shared<PsfexPsf>(*this);
92 }
93 
95 PsfexPsf::resized(int width, int height) const {
96  throw LSST_EXCEPT(pex::exceptions::LogicError, "Not Implemented");
97 }
98 
100 PsfexPsf::getKernel(afw::geom::Point2D position) const
101 {
102  double pos[MAXCONTEXT];
103  int const ndim = _context.size();
104  if (ndim != 2) { // we're only handling spatial variation for now
106  str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
107  % ndim));
108 
109  }
110  // where we want to evaluate the basis function's weights
111  if (!std::isfinite(position[0])) {
112  position = _averagePosition;
113  }
114 
115  for (int i = 0; i < ndim; ++i) {
116  pos[i] = (position[i] - _context[i].first)/_context[i].second;
117  }
118 
119  poly_func(_poly, pos); // evaluate polynomial
120 
121  int const w = _size[0], h = _size[1];
122  std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
123  /*
124  * Create a fixed Kernel out of each component, and then create a LinearCombinationKernel from them
125  */
126  const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
127  afw::math::KernelList kernels; kernels.reserve(nbasis); // the insides of the LinearCombinationKernel
128  std::vector<double> weights; weights.reserve(nbasis);
129 
130  float const vigstep = 1/_pixstep;
131  float const dx = 0.0, dy = 0.0;
132 
133  afw::geom::Box2I bbox = _doComputeBBox(position, afw::geom::Point2D(0, 0));
134  afw::detection::Psf::Image kim(bbox); // a basis function image, to be copied into a FixedKernel
135 
136  int sampleW = bbox.getWidth();
137  int sampleH = bbox.getHeight();
138 
139  std::vector<float> sampledBasis(sampleW*sampleH);
140 
141  for (int i = 0; i != nbasis; ++i) {
142  /*
143  * Resample the basis function onto the output resolution (and potentially subpixel offset)
144  */
145  vignet_resample(const_cast<float *>(&_comp[i*w*h]), w, h,
146  &sampledBasis[0], sampleW, sampleH,
147  -dx*vigstep, -dy*vigstep, vigstep, 1.0);
148  //
149  // And copy it into place
150  //
151  {
152  float *pl = &sampledBasis[0];
153  for (int y = 0; y != sampleH; ++y) {
154  for (int x = 0; x != sampleW; ++x) {
155  kim(x, y) = *pl++;
156  }
157  }
158  }
159 
160  kernels.push_back(std::make_shared<afw::math::FixedKernel>(kim));
161  weights.push_back(_poly->basis[i]);
162  }
163 
164  _kernel = std::make_shared<afw::math::LinearCombinationKernel>(kernels, weights);
165 
166  return _kernel;
167 }
168 
170 PsfexPsf::doComputeImage(afw::geom::Point2D const & position,
171  afw::image::Color const & color) const {
172  return _doComputeImage(position, color, position);
173 }
174 
176 PsfexPsf::doComputeKernelImage(afw::geom::Point2D const& position,
177  afw::image::Color const& color) const
178 {
179  return _doComputeImage(position, color, afw::geom::Point2D(0, 0));
180 }
181 
182 afw::geom::Box2I PsfexPsf::doComputeBBox(afw::geom::Point2D const & position,
183  afw::image::Color const & color) const {
184  return _doComputeBBox(position, afw::geom::Point2D(0, 0));
185 }
186 
187 afw::geom::Box2I PsfexPsf::_doComputeBBox(afw::geom::Point2D const & position,
188  afw::geom::Point2D const & center) const {
189  int const w = _size[0], h = _size[1];
190  int sampleW = static_cast<int>(w*_pixstep);
191  int sampleH = static_cast<int>(h*_pixstep);
192 
193  // Ensure that sizes are odd
194  if (sampleW % 2 == 0) sampleW += 1;
195  if (sampleH % 2 == 0) sampleH += 1;
196 
197  float dx = center[0] - static_cast<int>(center[0]);
198  float dy = center[1] - static_cast<int>(center[1]);
199 
200  if (dx > 0.5) dx -= 1.0;
201  if (dy > 0.5) dy -= 1.0;
202  // N.b. center[0] - dx == (int)center[x] until we reduced dx to (-0.5, 0.5].
203  // The + 0.5 is to handle floating point imprecision in this calculation
204  afw::geom::Box2I bbox(afw::geom::Point2I(static_cast<int>(center[0] - dx + 0.5) - sampleW/2,
205  static_cast<int>(center[1] - dy + 0.5) - sampleH/2),
206  afw::geom::Extent2I(sampleW, sampleH));
207  return bbox;
208 }
209 
211 PsfexPsf::_doComputeImage(afw::geom::Point2D const& position,
212  afw::image::Color const& color,
213  afw::geom::Point2D const& center
214  ) const
215 {
216  double pos[MAXCONTEXT];
217  int const ndim = _context.size();
218  if (ndim != 2) { // we're only handling spatial variation for now
220  str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
221  % ndim));
222 
223  }
224 
225  for (int i = 0; i < ndim; ++i) {
226  pos[i] = (position[i] - _context[i].first)/_context[i].second;
227  }
228 
229  poly_func(_poly, pos); // evaluate polynomial
230 
231  int const w = _size[0], h = _size[1];
232  std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
233  const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
234 
235  /* Sum each component */
236  int const npix = w*h;
237  for (int i = 0; i != nbasis; ++i) {
238  float *pl = &fullresIm[0];
239  float const fac = _poly->basis[i];
240  float const *ppc = &_comp[i*w*h];
241 
242  for (int j = 0; j != npix; ++j) {
243  pl[j] += fac*ppc[j];
244  }
245  }
246  /*
247  * We now have the image reconstructed at internal resolution; resample it onto the output resolution
248  * and subpixel offset
249  */
250  float const vigstep = 1/_pixstep;
251  float dx = center[0] - static_cast<int>(center[0]);
252  float dy = center[1] - static_cast<int>(center[1]);
253  if (dx > 0.5) dx -= 1.0;
254  if (dy > 0.5) dy -= 1.0;
255  //
256  // And copy it into place
257  //
258  afw::geom::Box2I bbox = _doComputeBBox(position, center);
259  PTR(afw::detection::Psf::Image) im = std::make_shared<afw::detection::Psf::Image>(bbox);
260 
261  int sampleW = bbox.getWidth();
262  int sampleH = bbox.getHeight();
263 
264  std::vector<float> sampledIm(sampleW*sampleH);
265 
266  vignet_resample(&fullresIm[0], w, h,
267  &sampledIm[0], sampleW, sampleH,
268  -dx*vigstep, -dy*vigstep, vigstep, 1.0);
269  {
270  float *pl = &sampledIm[0];
271  float const sum = std::accumulate(pl, pl + sampleW*sampleH, static_cast<float>(0));
272  for (int y = 0; y != sampleH; ++y) {
273  for (int x = 0; x != sampleW; ++x) {
274  (*im)(x, y) = *pl++/sum;
275  }
276  }
277  }
278 
279  return im;
280 }
281 
282 /************************************************************************************************************/
283 /*
284  * All the rest of this file handles persistence to FITS files
285  */
286 namespace table = afw::table;
287 
288 namespace {
289 
290 class PsfexPsfSchema1 {
291 public:
292  PsfexPsfSchema1() :
293  schema(),
294  ndim(schema.addField<int>("ndim", "Number of elements in group")),
295  ngroup(schema.addField<int>("ngroup", "Number of elements in degree")),
296  ncoeff(schema.addField<int>("ncoeff", "Number of coefficients")),
297 
298  _size_size(schema.addField<int>("_size_size", "Size of _size array")),
299  _comp_size(schema.addField<int>("_comp_size", "Size of _comp array")),
300  _context_size(schema.addField<int>("_context_size", "Size of _context array")),
301  // Other scalars
302  averagePosition(afw::table::PointKey<double>::addFields(schema,"averagePosition","average position of stars used to make the PSF","pixel")),
303  _pixstep(schema.addField<float>("_pixstep", "oversampling", "pixel"))
304  {
305  ;
306  }
307 
308  table::Schema schema;
309 
310  // Sizes in _poly
311  table::Key<int> ndim;
312  table::Key<int> ngroup;
313  table::Key<int> ncoeff;
314  // Sizes of vectors
315  table::Key<int> _size_size;
316  table::Key<int> _comp_size;
317  table::Key<int> _context_size;
318  // Other scalars
319  table::PointKey<double> averagePosition;
320  table::Key<float> _pixstep;
321 };
322 
323 class PsfexPsfSchema2 {
324 public:
325  PsfexPsfSchema2(int const ndim, int const ngroup, int const ncoeff,
326  int size_size, int comp_size, int context_size) :
327 
328  schema(),
329  group(schema.addField<table::Array<int> >("group", "Groups (of coefficients?)", ndim)),
330  degree(schema.addField<table::Array<int> >("degree", "Degree in each group", ngroup)),
331  basis(schema.addField<table::Array<double> >("basis", "Values of the basis functions", ncoeff)),
332  coeff(schema.addField<table::Array<double> >("coeff", "Polynomial coefficients", ncoeff)),
333  _size(schema.addField<table::Array<int> >("_size", "PSF dimensions", size_size)),
334  _comp(schema.addField<table::Array<float> >("_comp", "Complete pixel data", comp_size)),
335  _context_first( schema.addField<table::Array<double> >("_context_first",
336  "Offset to apply to context data", context_size)),
337  _context_second(schema.addField<table::Array<double> >("_context_second",
338  "Scale to apply to context data", context_size))
339  {
340  ;
341  }
342 
343  table::Schema schema;
344  // _poly
345  table::Key<table::Array<int> > group; // len(group) == ndim
346  table::Key<table::Array<int> > degree; // len(degree) == ngroup
347  table::Key<table::Array<double> > basis; // len(basis) == ncoeff
348  table::Key<table::Array<double> > coeff; // len(coeff) == ncoeff
349  // vectors
350  table::Key<table::Array<int> > _size;
351  table::Key<table::Array<float> > _comp;
352  table::Key<table::Array<double> > _context_first;
353  table::Key<table::Array<double> > _context_second;
354 };
355 
356 std::string getPsfexPsfPersistenceName() { return "PsfexPsf"; }
357 
358 } // anonymous
359 
360 /************************************************************************************************************/
361 
362 namespace detail { // PsfexPsfFactory needs to be a friend of PsfexPsf
363 class PsfexPsfFactory : public table::io::PersistableFactory {
364 public:
365 
366 virtual PTR(table::io::Persistable)
367 read(InputArchive const & archive, CatalogVector const & catalogs) const {
368  LSST_ARCHIVE_ASSERT(catalogs.size() == 2u);
369 
370  PTR(PsfexPsf) result(new PsfexPsf());
371 
372  int ndim, ngroup, ncoeff;
373  int size_size, comp_size, context_size;
374  {
375  PsfexPsfSchema1 const & keys = PsfexPsfSchema1();
376  LSST_ARCHIVE_ASSERT(catalogs[0].size() == 1u);
377  LSST_ARCHIVE_ASSERT(catalogs[0].getSchema() == keys.schema);
378  table::BaseRecord const & record = catalogs[0].front();
379 
380  // fields in _poly
381  ndim = record.get(keys.ndim);
382  ngroup = record.get(keys.ngroup);
383  ncoeff = record.get(keys.ncoeff);
384  // Other scalars
385  result->_averagePosition = record.get(keys.averagePosition);
386  result->_pixstep = record.get(keys._pixstep);
387  // sizes of vectors
388  size_size = record.get(keys._size_size);
389  comp_size = record.get(keys._comp_size);
390  context_size = record.get(keys._context_size);
391  }
392  // Now we can read the data
393  {
394  PsfexPsfSchema2 const keys(ndim, ngroup, ncoeff,
395  size_size, comp_size, context_size);
396 
397  LSST_ARCHIVE_ASSERT(catalogs[1].size() == 1u);
398  LSST_ARCHIVE_ASSERT(catalogs[1].getSchema() == keys.schema);
399  table::BaseRecord const & record = catalogs[1].front();
400 
401  // _poly
403  {
404  int const *begin = record.getElement(keys.group);
405  group.assign(begin, begin + ndim);
406 
407  for (int i = 0; i != ndim; ++i) {
408  ++group[i]; // poly_init subtracts 1 from each element. Sigh.
409  }
410  }
412  {
413  int const *begin = record.getElement(keys.degree);
414  degree.assign(begin, begin + ngroup);
415  }
416  result->_poly = poly_init(&group[0], group.size(), &degree[0], degree.size());
417  LSST_ARCHIVE_ASSERT(result->_poly->ncoeff == ncoeff);
418 
419  {
420  double const *begin = record.getElement(keys.basis);
421  std::copy(begin, begin + ncoeff, result->_poly->basis);
422  }
423  {
424  double const *begin = record.getElement(keys.coeff);
425  std::copy(begin, begin + ncoeff, result->_poly->coeff);
426  }
427  // vectors
428  {
429  int const *begin = record.getElement(keys._size);
430  result->_size.assign(begin, begin + size_size);
431  }
432  {
433  float const *begin = record.getElement(keys._comp);
434  result->_comp.assign(begin, begin + comp_size);
435  }
436  {
437  double const *begin1 = record.getElement(keys._context_first);
438  double const *begin2 = record.getElement(keys._context_second);
439  result->_context.resize(context_size);
440  for (int i = 0; i != context_size; ++i) {
441  result->_context[i].first = begin1[i];
442  result->_context[i].second = begin2[i];
443  }
444  }
445  }
446 
447  return result;
448 }
449 
450 explicit PsfexPsfFactory(std::string const & name) : table::io::PersistableFactory(name) {}
451 
452 };
453 }
454 
455 namespace {
456  detail::PsfexPsfFactory registration(getPsfexPsfPersistenceName());
457 }
458 
459 /************************************************************************************************************/
460 
461 std::string PsfexPsf::getPythonModule() const { return "lsst.meas.extensions.psfex"; }
462 
463 std::string PsfexPsf::getPersistenceName() const { return getPsfexPsfPersistenceName(); }
464 
466  // First write the dimensions of the various arrays to an HDU, so we can construct the second
467  // HDU using them when we read the Psf back
468  {
469  PsfexPsfSchema1 const keys;
470  afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
471  PTR(afw::table::BaseRecord) record = cat.addNew();
472 
473  // Sizes in _poly
474  record->set(keys.ndim, _poly->ndim);
475  record->set(keys.ngroup, _poly->ngroup);
476  record->set(keys.ncoeff, _poly->ncoeff);
477  // Other scalars
478  record->set(keys.averagePosition, _averagePosition);
479  record->set(keys._pixstep, _pixstep);
480  record->set(keys._size_size, _size.size());
481  record->set(keys._comp_size, _comp.size());
482  record->set(keys._context_size, _context.size());
483 
484  handle.saveCatalog(cat);
485  }
486  // Now we can write the data
487  {
488  PsfexPsfSchema2 const keys(_poly->ndim, _poly->ngroup, _poly->ncoeff,
489  _size.size(), _comp.size(), _context.size());
490  afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
491  // _poly
492  PTR(afw::table::BaseRecord) record = cat.addNew();
493  {
494  int *begin = record->getElement(keys.group);
495  std::copy(_poly->group, _poly->group + _poly->ndim, begin);
496  }
497  {
498  int *begin = record->getElement(keys.degree);
499  std::copy(_poly->degree, _poly->degree + _poly->ngroup, begin);
500  }
501  {
502  double *begin = record->getElement(keys.basis);
503  std::copy(_poly->basis, _poly->basis + _poly->ncoeff, begin);
504  }
505  {
506  double *begin = record->getElement(keys.coeff);
507  std::copy(_poly->coeff, _poly->coeff + _poly->ncoeff, begin);
508  }
509  // vectors
510  {
511  int *begin = record->getElement(keys._size);
512  std::copy(_size.begin(), _size.end(), begin);
513  }
514  {
515  float *begin = record->getElement(keys._comp);
516  std::copy(_comp.begin(), _comp.end(), begin);
517  }
518  {
519  double *begin1 = record->getElement(keys._context_first);
520  double *begin2 = record->getElement(keys._context_second);
521  for (unsigned int i = 0; i != _context.size(); ++i) {
522  begin1[i] = _context[i].first;
523  begin2[i] = _context[i].second;
524  }
525  }
526 
527  handle.saveCatalog(cat);
528  }
529 }
530 
531 }}}}
table::Key< int > ncoeff
Definition: PsfexPsf.cc:313
T copy(T... args)
#define LSST_ARCHIVE_ASSERT(EXPR)
table::Key< int > ndim
Definition: PsfexPsf.cc:311
#define PTR(...)
table::Key< table::Array< int > > degree
Definition: PsfexPsf.cc:346
Field< T >::Value get(Key< T > const &key) const
table::Key< int > _size_size
Definition: PsfexPsf.cc:315
virtual boost::shared_ptr< lsst::afw::detection::Psf > clone() const
Polymorphic deep copy; should usually be unnecessary as Psfs are immutable.x.
Definition: PsfexPsf.cc:90
STL class.
T push_back(T... args)
Represent a PSF as a linear combination of PSFEX (== Karhunen-Loeve) basis functions.
Definition: PsfexPsf.h:39
virtual boost::shared_ptr< afw::detection::Psf > resized(int width, int height) const
Return a clone with specified kernel dimensions.
Definition: PsfexPsf.cc:95
double w
table::Key< table::Array< double > > _context_second
Definition: PsfexPsf.cc:353
ImagePsf(bool isFixed=false)
def keys(self)
T isfinite(T... args)
double x
table::Key< table::Array< int > > _size
Definition: PsfexPsf.cc:350
table::Key< table::Array< double > > coeff
Definition: PsfexPsf.cc:348
std::unique_ptr< SchemaItem< U > > result
BaseCatalog makeCatalog(Schema const &schema)
Field< T >::Element * getElement(Key< T > const &key)
table::Key< int > _comp_size
Definition: PsfexPsf.cc:316
second
T size(T... args)
boost::shared_ptr< lsst::afw::math::LinearCombinationKernel const > getKernel(lsst::afw::geom::Point2D=lsst::afw::geom::Point2D(std::numeric_limits< double >::quiet_NaN())) const
Return the PSF&#39;s basis functions as a spatially-invariant LinearCombinationKernel with unit weights...
Definition: PsfexPsf.cc:100
#define LSST_EXCEPT(type,...)
T assign(T... args)
table::Key< table::Array< double > > basis
Definition: PsfexPsf.cc:347
table::Key< table::Array< float > > _comp
Definition: PsfexPsf.cc:351
void write(lsst::afw::table::io::OutputArchiveHandle &handle) const
Definition: PsfexPsf.cc:465
table::Key< table::Array< double > > _context_first
Definition: PsfexPsf.cc:352
table::Schema schema
Definition: PsfexPsf.cc:308
table::Box2IKey bbox
table::Key< int > _context_size
Definition: PsfexPsf.cc:317
void saveCatalog(BaseCatalog const &catalog)
table::Key< table::Array< int > > group
Definition: PsfexPsf.cc:345
T accumulate(T... args)
table::Key< float > _pixstep
Definition: PsfexPsf.cc:320
table::PointKey< double > averagePosition
Definition: PsfexPsf.cc:319
int y
T reserve(T... args)
table::Key< int > ngroup
Definition: PsfexPsf.cc:312
std::shared_ptr< BaseRecord > addNew()