LSST Applications  21.0.0-147-g0e635eb1+1acddb5be5,22.0.0+052faf71bd,22.0.0+1ea9a8b2b2,22.0.0+6312710a6c,22.0.0+729191ecac,22.0.0+7589c3a021,22.0.0+9f079a9461,22.0.1-1-g7d6de66+b8044ec9de,22.0.1-1-g87000a6+536b1ee016,22.0.1-1-g8e32f31+6312710a6c,22.0.1-10-gd060f87+016f7cdc03,22.0.1-12-g9c3108e+df145f6f68,22.0.1-16-g314fa6d+c825727ab8,22.0.1-19-g93a5c75+d23f2fb6d8,22.0.1-19-gb93eaa13+aab3ef7709,22.0.1-2-g8ef0a89+b8044ec9de,22.0.1-2-g92698f7+9f079a9461,22.0.1-2-ga9b0f51+052faf71bd,22.0.1-2-gac51dbf+052faf71bd,22.0.1-2-gb66926d+6312710a6c,22.0.1-2-gcb770ba+09e3807989,22.0.1-20-g32debb5+b8044ec9de,22.0.1-23-gc2439a9a+fb0756638e,22.0.1-3-g496fd5d+09117f784f,22.0.1-3-g59f966b+1e6ba2c031,22.0.1-3-g849a1b8+f8b568069f,22.0.1-3-gaaec9c0+c5c846a8b1,22.0.1-32-g5ddfab5d3+60ce4897b0,22.0.1-4-g037fbe1+64e601228d,22.0.1-4-g8623105+b8044ec9de,22.0.1-5-g096abc9+d18c45d440,22.0.1-5-g15c806e+57f5c03693,22.0.1-7-gba73697+57f5c03693,master-g6e05de7fdc+c1283a92b8,master-g72cdda8301+729191ecac,w.2021.39
LSST Data Management Base Package
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"
40 extern "C" {
41 #include "vignet.h"
42 }
46 
47 namespace lsst {
48 namespace afw {
49 namespace table {
50 namespace io {
51 
54 
55 } // namespace io
56 } // namespace table
57 } // namespace afw
58 namespace meas {
59 namespace extensions {
60 namespace psfex {
61 
62 namespace afw = lsst::afw;
63 
64 PsfexPsf::PsfexPsf(
65  lsst::meas::extensions::psfex::Psf const& psf,
67  ) : ImagePsf(), _averagePosition(averagePosition),
68  _size(psf.impl->dim),
69  _comp(psf.impl->npix),
70  _context(psf.impl->poly->ndim)
71 
72 {
73  _poly = poly_copy(psf.impl->poly);
74 
75  _pixstep = psf.impl->pixstep;
76 
77  std::copy(psf.impl->size, psf.impl->size + psf.impl->dim, _size.begin());
78 
79  std::copy(psf.impl->comp, psf.impl->comp + psf.impl->npix, _comp.begin());
80 
81  for (int i = 0; i != psf.impl->poly->ndim; ++i) {
82  _context[i].first = psf.impl->contextoffset[i];
83  _context[i].second = psf.impl->contextscale[i];
84  }
85 }
86 
87  PsfexPsf::PsfexPsf() : ImagePsf(),
88  _averagePosition(geom::Point2I(0, 0)),
89  _poly(0),
90  _pixstep(0.0),
91  _size(),
92  _comp(),
93  _context()
94 {
95  ;
96 }
97 
99 {
100  poly_end(_poly);
101 }
102 
105  return std::make_shared<PsfexPsf>(*this);
106 }
107 
109 PsfexPsf::resized(int width, int height) const {
110  throw LSST_EXCEPT(pex::exceptions::LogicError, "Not Implemented");
111 }
112 
115 {
116  double pos[MAXCONTEXT];
117  int const ndim = _context.size();
118  if (ndim != 2) { // we're only handling spatial variation for now
120  str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
121  % ndim));
122 
123  }
124  // where we want to evaluate the basis function's weights
125  if (!std::isfinite(position[0])) {
126  position = _averagePosition;
127  }
128 
129  for (int i = 0; i < ndim; ++i) {
130  pos[i] = (position[i] - _context[i].first)/_context[i].second;
131  }
132 
133  poly_func(_poly, pos); // evaluate polynomial
134 
135  int const w = _size[0], h = _size[1];
136  std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
137  /*
138  * Create a fixed Kernel out of each component, and then create a LinearCombinationKernel from them
139  */
140  const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
141  afw::math::KernelList kernels; kernels.reserve(nbasis); // the insides of the LinearCombinationKernel
142  std::vector<double> weights; weights.reserve(nbasis);
143 
144  float const vigstep = 1/_pixstep;
145  float const dx = 0.0, dy = 0.0;
146 
147  geom::Box2I bbox = _doComputeBBox(position, geom::Point2D(0, 0));
148  afw::detection::Psf::Image kim(bbox); // a basis function image, to be copied into a FixedKernel
149 
150  int sampleW = bbox.getWidth();
151  int sampleH = bbox.getHeight();
152 
153  std::vector<float> sampledBasis(sampleW*sampleH);
154 
155  for (int i = 0; i != nbasis; ++i) {
156  /*
157  * Resample the basis function onto the output resolution (and potentially subpixel offset)
158  */
159  vignet_resample(const_cast<float *>(&_comp[i*w*h]), w, h,
160  &sampledBasis[0], sampleW, sampleH,
161  -dx*vigstep, -dy*vigstep, vigstep, 1.0);
162  //
163  // And copy it into place
164  //
165  {
166  float *pl = &sampledBasis[0];
167  for (int y = 0; y != sampleH; ++y) {
168  for (int x = 0; x != sampleW; ++x) {
169  kim(x, y) = *pl++;
170  }
171  }
172  }
173 
174  kernels.push_back(std::make_shared<afw::math::FixedKernel>(kim));
175  weights.push_back(_poly->basis[i]);
176  }
177 
178  _kernel = std::make_shared<afw::math::LinearCombinationKernel>(kernels, weights);
179 
180  return _kernel;
181 }
182 
184 PsfexPsf::doComputeImage(geom::Point2D const & position,
185  afw::image::Color const & color) const {
186  return _doComputeImage(position, color, position);
187 }
188 
190 PsfexPsf::doComputeKernelImage(geom::Point2D const& position,
191  afw::image::Color const& color) const
192 {
193  return _doComputeImage(position, color, geom::Point2D(0, 0));
194 }
195 
196 geom::Box2I PsfexPsf::doComputeBBox(geom::Point2D const & position,
197  afw::image::Color const & color) const {
198  return _doComputeBBox(position, geom::Point2D(0, 0));
199 }
200 
201 geom::Box2I PsfexPsf::_doComputeBBox(geom::Point2D const & position,
202  geom::Point2D const & center) const {
203  int const w = _size[0], h = _size[1];
204  int sampleW = static_cast<int>(w*_pixstep);
205  int sampleH = static_cast<int>(h*_pixstep);
206 
207  // Ensure that sizes are odd
208  if (sampleW % 2 == 0) sampleW += 1;
209  if (sampleH % 2 == 0) sampleH += 1;
210 
211  float dx = center[0] - static_cast<int>(center[0]);
212  float dy = center[1] - static_cast<int>(center[1]);
213 
214  if (dx > 0.5) dx -= 1.0;
215  if (dy > 0.5) dy -= 1.0;
216  // N.b. center[0] - dx == (int)center[x] until we reduced dx to (-0.5, 0.5].
217  // The + 0.5 is to handle floating point imprecision in this calculation
218  geom::Box2I bbox(geom::Point2I(static_cast<int>(center[0] - dx + 0.5) - sampleW/2,
219  static_cast<int>(center[1] - dy + 0.5) - sampleH/2),
220  geom::Extent2I(sampleW, sampleH));
221  return bbox;
222 }
223 
226  afw::image::Color const& color,
227  geom::Point2D const& center
228  ) const
229 {
230  double pos[MAXCONTEXT];
231  int const ndim = _context.size();
232  if (ndim != 2) { // we're only handling spatial variation for now
234  str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
235  % ndim));
236 
237  }
238 
239  for (int i = 0; i < ndim; ++i) {
240  pos[i] = (position[i] - _context[i].first)/_context[i].second;
241  }
242 
243  poly_func(_poly, pos); // evaluate polynomial
244 
245  int const w = _size[0], h = _size[1];
246  std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
247  const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
248 
249  /* Sum each component */
250  int const npix = w*h;
251  for (int i = 0; i != nbasis; ++i) {
252  float *pl = &fullresIm[0];
253  float const fac = _poly->basis[i];
254  float const *ppc = &_comp[i*w*h];
255 
256  for (int j = 0; j != npix; ++j) {
257  pl[j] += fac*ppc[j];
258  }
259  }
260  /*
261  * We now have the image reconstructed at internal resolution; resample it onto the output resolution
262  * and subpixel offset
263  */
264  float const vigstep = 1/_pixstep;
265  float dx = center[0] - static_cast<int>(center[0]);
266  float dy = center[1] - static_cast<int>(center[1]);
267  if (dx > 0.5) dx -= 1.0;
268  if (dy > 0.5) dy -= 1.0;
269  //
270  // And copy it into place
271  //
272  geom::Box2I bbox = _doComputeBBox(position, center);
273  std::shared_ptr<afw::detection::Psf::Image> im = std::make_shared<afw::detection::Psf::Image>(bbox);
274 
275  int sampleW = bbox.getWidth();
276  int sampleH = bbox.getHeight();
277 
278  std::vector<float> sampledIm(sampleW*sampleH);
279 
280  vignet_resample(&fullresIm[0], w, h,
281  &sampledIm[0], sampleW, sampleH,
282  -dx*vigstep, -dy*vigstep, vigstep, 1.0);
283  {
284  float *pl = &sampledIm[0];
285  float const sum = std::accumulate(pl, pl + sampleW*sampleH, static_cast<float>(0));
286  for (int y = 0; y != sampleH; ++y) {
287  for (int x = 0; x != sampleW; ++x) {
288  (*im)(x, y) = *pl++/sum;
289  }
290  }
291  }
292 
293  return im;
294 }
295 
296 /************************************************************************************************************/
297 /*
298  * All the rest of this file handles persistence to FITS files
299  */
300 namespace table = afw::table;
301 
302 namespace {
303 
304 class PsfexPsfSchema1 {
305 public:
306  PsfexPsfSchema1() :
307  schema(),
308  ndim(schema.addField<int>("ndim", "Number of elements in group")),
309  ngroup(schema.addField<int>("ngroup", "Number of elements in degree")),
310  ncoeff(schema.addField<int>("ncoeff", "Number of coefficients")),
311 
312  _size_size(schema.addField<int>("_size_size", "Size of _size array")),
313  _comp_size(schema.addField<int>("_comp_size", "Size of _comp array")),
314  _context_size(schema.addField<int>("_context_size", "Size of _context array")),
315  // Other scalars
316  averagePosition(afw::table::PointKey<double>::addFields(schema,"averagePosition","average position of stars used to make the PSF","pixel")),
317  _pixstep(schema.addField<float>("_pixstep", "oversampling", "pixel"))
318  {
319  ;
320  }
321 
322  table::Schema schema;
323 
324  // Sizes in _poly
325  table::Key<int> ndim;
326  table::Key<int> ngroup;
327  table::Key<int> ncoeff;
328  // Sizes of vectors
329  table::Key<int> _size_size;
330  table::Key<int> _comp_size;
331  table::Key<int> _context_size;
332  // Other scalars
333  table::PointKey<double> averagePosition;
334  table::Key<float> _pixstep;
335 };
336 
337 class PsfexPsfSchema2 {
338 public:
339  PsfexPsfSchema2(int const ndim, int const ngroup, int const ncoeff,
340  int size_size, int comp_size, int context_size) :
341 
342  schema(),
343  group(schema.addField<table::Array<int> >("group", "Groups (of coefficients?)", ndim)),
344  degree(schema.addField<table::Array<int> >("degree", "Degree in each group", ngroup)),
345  basis(schema.addField<table::Array<double> >("basis", "Values of the basis functions", ncoeff)),
346  coeff(schema.addField<table::Array<double> >("coeff", "Polynomial coefficients", ncoeff)),
347  _size(schema.addField<table::Array<int> >("_size", "PSF dimensions", size_size)),
348  _comp(schema.addField<table::Array<float> >("_comp", "Complete pixel data", comp_size)),
349  _context_first( schema.addField<table::Array<double> >("_context_first",
350  "Offset to apply to context data", context_size)),
351  _context_second(schema.addField<table::Array<double> >("_context_second",
352  "Scale to apply to context data", context_size))
353  {
354  ;
355  }
356 
357  table::Schema schema;
358  // _poly
359  table::Key<table::Array<int> > group; // len(group) == ndim
360  table::Key<table::Array<int> > degree; // len(degree) == ngroup
361  table::Key<table::Array<double> > basis; // len(basis) == ncoeff
362  table::Key<table::Array<double> > coeff; // len(coeff) == ncoeff
363  // vectors
364  table::Key<table::Array<int> > _size;
365  table::Key<table::Array<float> > _comp;
366  table::Key<table::Array<double> > _context_first;
367  table::Key<table::Array<double> > _context_second;
368 };
369 
370 std::string getPsfexPsfPersistenceName() { return "PsfexPsf"; }
371 
372 } // anonymous
373 
374 /************************************************************************************************************/
375 
376 namespace detail { // PsfexPsfFactory needs to be a friend of PsfexPsf
377 class PsfexPsfFactory : public table::io::PersistableFactory {
378 public:
379 
381 read(InputArchive const & archive, CatalogVector const & catalogs) const {
382  LSST_ARCHIVE_ASSERT(catalogs.size() == 2u);
383 
385 
386  int ndim, ngroup, ncoeff;
387  int size_size, comp_size, context_size;
388  {
389  PsfexPsfSchema1 const & keys = PsfexPsfSchema1();
390  LSST_ARCHIVE_ASSERT(catalogs[0].size() == 1u);
391  LSST_ARCHIVE_ASSERT(catalogs[0].getSchema() == keys.schema);
392  table::BaseRecord const & record = catalogs[0].front();
393 
394  // fields in _poly
395  ndim = record.get(keys.ndim);
396  ngroup = record.get(keys.ngroup);
397  ncoeff = record.get(keys.ncoeff);
398  // Other scalars
399  result->_averagePosition = record.get(keys.averagePosition);
400  result->_pixstep = record.get(keys._pixstep);
401  // sizes of vectors
402  size_size = record.get(keys._size_size);
403  comp_size = record.get(keys._comp_size);
404  context_size = record.get(keys._context_size);
405  }
406  // Now we can read the data
407  {
408  PsfexPsfSchema2 const keys(ndim, ngroup, ncoeff,
409  size_size, comp_size, context_size);
410 
411  LSST_ARCHIVE_ASSERT(catalogs[1].size() == 1u);
412  LSST_ARCHIVE_ASSERT(catalogs[1].getSchema() == keys.schema);
413  table::BaseRecord const & record = catalogs[1].front();
414 
415  // _poly
417  {
418  int const *begin = record.getElement(keys.group);
419  group.assign(begin, begin + ndim);
420 
421  for (int i = 0; i != ndim; ++i) {
422  ++group[i]; // poly_init subtracts 1 from each element. Sigh.
423  }
424  }
426  {
427  int const *begin = record.getElement(keys.degree);
428  degree.assign(begin, begin + ngroup);
429  }
430  result->_poly = poly_init(&group[0], group.size(), &degree[0], degree.size());
431  LSST_ARCHIVE_ASSERT(result->_poly->ncoeff == ncoeff);
432 
433  {
434  double const *begin = record.getElement(keys.basis);
435  std::copy(begin, begin + ncoeff, result->_poly->basis);
436  }
437  {
438  double const *begin = record.getElement(keys.coeff);
439  std::copy(begin, begin + ncoeff, result->_poly->coeff);
440  }
441  // vectors
442  {
443  int const *begin = record.getElement(keys._size);
444  result->_size.assign(begin, begin + size_size);
445  }
446  {
447  float const *begin = record.getElement(keys._comp);
448  result->_comp.assign(begin, begin + comp_size);
449  }
450  {
451  double const *begin1 = record.getElement(keys._context_first);
452  double const *begin2 = record.getElement(keys._context_second);
453  result->_context.resize(context_size);
454  for (int i = 0; i != context_size; ++i) {
455  result->_context[i].first = begin1[i];
456  result->_context[i].second = begin2[i];
457  }
458  }
459  }
460 
461  return result;
462 }
463 
464 explicit PsfexPsfFactory(std::string const & name) : table::io::PersistableFactory(name) {}
465 
466 };
467 }
468 
469 namespace {
470  detail::PsfexPsfFactory registration(getPsfexPsfPersistenceName());
471 }
472 
473 /************************************************************************************************************/
474 
475 std::string PsfexPsf::getPythonModule() const { return "lsst.meas.extensions.psfex"; }
476 
477 std::string PsfexPsf::getPersistenceName() const { return getPsfexPsfPersistenceName(); }
478 
480  // First write the dimensions of the various arrays to an HDU, so we can construct the second
481  // HDU using them when we read the Psf back
482  {
483  PsfexPsfSchema1 const keys;
484  afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
486 
487  // Sizes in _poly
488  record->set(keys.ndim, _poly->ndim);
489  record->set(keys.ngroup, _poly->ngroup);
490  record->set(keys.ncoeff, _poly->ncoeff);
491  // Other scalars
492  record->set(keys.averagePosition, _averagePosition);
493  record->set(keys._pixstep, _pixstep);
494  record->set(keys._size_size, _size.size());
495  record->set(keys._comp_size, _comp.size());
496  record->set(keys._context_size, _context.size());
497 
498  handle.saveCatalog(cat);
499  }
500  // Now we can write the data
501  {
502  PsfexPsfSchema2 const keys(_poly->ndim, _poly->ngroup, _poly->ncoeff,
503  _size.size(), _comp.size(), _context.size());
504  afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
505  // _poly
507  {
508  int *begin = record->getElement(keys.group);
509  std::copy(_poly->group, _poly->group + _poly->ndim, begin);
510  }
511  {
512  int *begin = record->getElement(keys.degree);
513  std::copy(_poly->degree, _poly->degree + _poly->ngroup, begin);
514  }
515  {
516  double *begin = record->getElement(keys.basis);
517  std::copy(_poly->basis, _poly->basis + _poly->ncoeff, begin);
518  }
519  {
520  double *begin = record->getElement(keys.coeff);
521  std::copy(_poly->coeff, _poly->coeff + _poly->ncoeff, begin);
522  }
523  // vectors
524  {
525  int *begin = record->getElement(keys._size);
526  std::copy(_size.begin(), _size.end(), begin);
527  }
528  {
529  float *begin = record->getElement(keys._comp);
530  std::copy(_comp.begin(), _comp.end(), begin);
531  }
532  {
533  double *begin1 = record->getElement(keys._context_first);
534  double *begin2 = record->getElement(keys._context_second);
535  for (unsigned int i = 0; i != _context.size(); ++i) {
536  begin1[i] = _context[i].first;
537  begin2[i] = _context[i].second;
538  }
539  }
540 
541  handle.saveCatalog(cat);
542  }
543 }
544 
545 }}}}
py::object result
Definition: _schema.cc:429
table::Key< std::string > name
Definition: Amplifier.cc:116
AmpInfoBoxKey bbox
Definition: Amplifier.cc:117
double x
#define LSST_EXCEPT(type,...)
Create an exception with a given type.
Definition: Exception.h:48
Utilities for persisting KernelPsf and subclasses thereof.
int y
Definition: SpanSet.cc:48
T accumulate(T... args)
#define LSST_ARCHIVE_ASSERT(EXPR)
An assertion macro used to validate the structure of an InputArchive.
Definition: Persistable.h:48
Basic LSST definitions.
T begin(T... args)
Describe the colour of a source.
Definition: Color.h:26
Base class for all records.
Definition: BaseRecord.h:31
Field< T >::Element * getElement(Key< T > const &key)
Return a pointer to the underlying elements of a field (non-const).
Definition: BaseRecord.h:93
Field< T >::Value get(Key< T > const &key) const
Return the value of a field for the given key.
Definition: BaseRecord.h:151
std::shared_ptr< RecordT > addNew()
Create a new record, add it to the end of the catalog, and return a pointer to it.
Definition: Catalog.h:490
An object passed to Persistable::write to allow it to persist itself.
void saveCatalog(BaseCatalog const &catalog)
Save a catalog in the archive.
BaseCatalog makeCatalog(Schema const &schema)
Return a new, empty catalog with the given schema.
static std::shared_ptr< T > dynamicCast(std::shared_ptr< Persistable > const &ptr)
Dynamically cast a shared_ptr.
Definition: Persistable.cc:18
An integer coordinate rectangle.
Definition: Box.h:55
Represent a PSF as a linear combination of PSFEX (== Karhunen-Loeve) basis functions.
Definition: PsfexPsf.h:40
std::shared_ptr< lsst::afw::math::LinearCombinationKernel const > getKernel(lsst::geom::Point2D=lsst::geom::Point2D(std::numeric_limits< double >::quiet_NaN())) const
Return the PSF's basis functions as a spatially-invariant LinearCombinationKernel with unit weights.
Definition: PsfexPsf.cc:114
void write(lsst::afw::table::io::OutputArchiveHandle &handle) const override
Write the object to one or more catalogs.
Definition: PsfexPsf.cc:479
std::shared_ptr< lsst::afw::detection::Psf > clone() const override
Polymorphic deep copy; should usually be unnecessary as Psfs are immutable.x.
Definition: PsfexPsf.cc:104
virtual std::shared_ptr< lsst::afw::detection::Psf::Image > _doComputeImage(lsst::geom::Point2D const &position, lsst::afw::image::Color const &color, lsst::geom::Point2D const &center) const
Compute an image of the Psf at the specified position/colour, at pixel position in the output image.
Definition: PsfexPsf.cc:225
std::shared_ptr< afw::detection::Psf > resized(int width, int height) const override
Return a clone with specified kernel dimensions.
Definition: PsfexPsf.cc:109
virtual std::shared_ptr< table::io::Persistable > read(InputArchive const &archive, CatalogVector const &catalogs) const
Definition: PsfexPsf.cc:381
Reports invalid arguments.
Definition: Runtime.h:66
Reports errors in the logical structure of the program.
Definition: Runtime.h:46
T copy(T... args)
T end(T... args)
T isfinite(T... args)
Low-level polynomials (including special polynomials) in C++.
Definition: Basis1d.h:26
Point< int, 2 > Point2I
Definition: Point.h:321
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
A base class for image defects.
T push_back(T... args)
T reserve(T... args)
T size(T... args)
double w
Definition: CoaddPsf.cc:69
table::Key< int > ncoeff
Definition: PsfexPsf.cc:327
table::Key< int > _context_size
Definition: PsfexPsf.cc:331
table::Key< int > ndim
Definition: PsfexPsf.cc:325
table::Key< int > ngroup
Definition: PsfexPsf.cc:326
table::Schema schema
Definition: PsfexPsf.cc:322
table::Key< table::Array< int > > _size
Definition: PsfexPsf.cc:364
table::Key< table::Array< int > > degree
Definition: PsfexPsf.cc:360
table::PointKey< double > averagePosition
Definition: PsfexPsf.cc:333
table::Key< table::Array< float > > _comp
Definition: PsfexPsf.cc:365
table::Key< float > _pixstep
Definition: PsfexPsf.cc:334
table::Key< table::Array< double > > _context_first
Definition: PsfexPsf.cc:366
table::Key< table::Array< double > > coeff
Definition: PsfexPsf.cc:362
table::Key< table::Array< double > > basis
Definition: PsfexPsf.cc:361
table::Key< table::Array< double > > _context_second
Definition: PsfexPsf.cc:367
table::Key< table::Array< int > > group
Definition: PsfexPsf.cc:359
table::Key< int > _comp_size
Definition: PsfexPsf.cc:330
table::Key< int > _size_size
Definition: PsfexPsf.cc:329
Key< int > psf
Definition: Exposure.cc:65