LSST Applications g0f08755f38+82efc23009,g12f32b3c4e+e7bdf1200e,g1653933729+a8ce1bb630,g1a0ca8cf93+50eff2b06f,g28da252d5a+52db39f6a5,g2bbee38e9b+37c5a29d61,g2bc492864f+37c5a29d61,g2cdde0e794+c05ff076ad,g3156d2b45e+41e33cbcdc,g347aa1857d+37c5a29d61,g35bb328faa+a8ce1bb630,g3a166c0a6a+37c5a29d61,g3e281a1b8c+fb992f5633,g414038480c+7f03dfc1b0,g41af890bb2+11b950c980,g5fbc88fb19+17cd334064,g6b1c1869cb+12dd639c9a,g781aacb6e4+a8ce1bb630,g80478fca09+72e9651da0,g82479be7b0+04c31367b4,g858d7b2824+82efc23009,g9125e01d80+a8ce1bb630,g9726552aa6+8047e3811d,ga5288a1d22+e532dc0a0b,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+37c5a29d61,gcf0d15dbbd+2acd6d4d48,gd7358e8bfb+778a810b6e,gda3e153d99+82efc23009,gda6a2b7d83+2acd6d4d48,gdaeeff99f8+1711a396fd,ge2409df99d+6b12de1076,ge79ae78c31+37c5a29d61,gf0baf85859+d0a5978c5a,gf3967379c6+4954f8c433,gfb92a5be7c+82efc23009,gfec2e1e490+2aaed99252,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
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"
40extern "C" {
41#include "vignet.h"
42}
46
47namespace lsst {
48namespace afw {
49namespace table {
50namespace io {
51
54
55} // namespace io
56} // namespace table
57} // namespace afw
58namespace meas {
59namespace extensions {
60namespace psfex {
61
62namespace afw = lsst::afw;
63
64PsfexPsf::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
109PsfexPsf::resized(int width, int height) const {
110 throw LSST_EXCEPT(pex::exceptions::LogicError, "Not Implemented");
111}
112
113int PsfexPsf::getNdim() const {
114 return _context.size();
115}
116
119{
120 double pos[MAXCONTEXT];
121 int const ndim = _context.size();
122 if (ndim != 2) { // we're only handling spatial variation for now
124 str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
125 % ndim));
126
127 }
128 // where we want to evaluate the basis function's weights
129 if (!std::isfinite(position[0])) {
130 position = _averagePosition;
131 }
132
133 for (int i = 0; i < ndim; ++i) {
134 pos[i] = (position[i] - _context[i].first)/_context[i].second;
135 }
136
137 poly_func(_poly, pos); // evaluate polynomial
138
139 int const w = _size[0], h = _size[1];
140 std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
141 /*
142 * Create a fixed Kernel out of each component, and then create a LinearCombinationKernel from them
143 */
144 const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
145 afw::math::KernelList kernels; kernels.reserve(nbasis); // the insides of the LinearCombinationKernel
146 std::vector<double> weights; weights.reserve(nbasis);
147
148 float const vigstep = 1/_pixstep;
149 float const dx = 0.0, dy = 0.0;
150
151 geom::Box2I bbox = _doComputeBBox(position, geom::Point2D(0, 0));
152 afw::detection::Psf::Image kim(bbox); // a basis function image, to be copied into a FixedKernel
153
154 int sampleW = bbox.getWidth();
155 int sampleH = bbox.getHeight();
156
157 std::vector<float> sampledBasis(sampleW*sampleH);
158
159 for (int i = 0; i != nbasis; ++i) {
160 /*
161 * Resample the basis function onto the output resolution (and potentially subpixel offset)
162 */
163 vignet_resample(const_cast<float *>(&_comp[i*w*h]), w, h,
164 &sampledBasis[0], sampleW, sampleH,
165 -dx*vigstep, -dy*vigstep, vigstep, 1.0);
166 //
167 // And copy it into place
168 //
169 {
170 float *pl = &sampledBasis[0];
171 for (int y = 0; y != sampleH; ++y) {
172 for (int x = 0; x != sampleW; ++x) {
173 kim(x, y) = *pl++;
174 }
175 }
176 }
177
178 kernels.push_back(std::make_shared<afw::math::FixedKernel>(kim));
179 weights.push_back(_poly->basis[i]);
180 }
181
182 _kernel = std::make_shared<afw::math::LinearCombinationKernel>(kernels, weights);
183
184 return _kernel;
185}
186
189 afw::image::Color const & color) const {
190 return _doComputeImage(position, color, position);
191}
192
195 afw::image::Color const& color) const
196{
197 return _doComputeImage(position, color, geom::Point2D(0, 0));
198}
199
201 afw::image::Color const & color) const {
202 return _doComputeBBox(position, geom::Point2D(0, 0));
203}
204
205geom::Box2I PsfexPsf::_doComputeBBox(geom::Point2D const & position,
206 geom::Point2D const & center) const {
207 int const w = _size[0], h = _size[1];
208 int sampleW = static_cast<int>(w*_pixstep);
209 int sampleH = static_cast<int>(h*_pixstep);
210
211 // Ensure that sizes are odd
212 if (sampleW % 2 == 0) sampleW += 1;
213 if (sampleH % 2 == 0) sampleH += 1;
214
215 float dx = center[0] - static_cast<int>(center[0]);
216 float dy = center[1] - static_cast<int>(center[1]);
217
218 if (dx > 0.5) dx -= 1.0;
219 if (dy > 0.5) dy -= 1.0;
220 // N.b. center[0] - dx == (int)center[x] until we reduced dx to (-0.5, 0.5].
221 // The + 0.5 is to handle floating point imprecision in this calculation
222 geom::Box2I bbox(geom::Point2I(static_cast<int>(center[0] - dx + 0.5) - sampleW/2,
223 static_cast<int>(center[1] - dy + 0.5) - sampleH/2),
224 geom::Extent2I(sampleW, sampleH));
225 return bbox;
226}
227
230 afw::image::Color const& color,
231 geom::Point2D const& center
232 ) const
233{
234 double pos[MAXCONTEXT];
235 int const ndim = _context.size();
236 if (ndim != 2) { // we're only handling spatial variation for now
238 str(boost::format("Only spatial variation (ndim == 2) is supported; saw %d")
239 % ndim));
240
241 }
242
243 for (int i = 0; i < ndim; ++i) {
244 pos[i] = (position[i] - _context[i].first)/_context[i].second;
245 }
246
247 poly_func(_poly, pos); // evaluate polynomial
248
249 int const w = _size[0], h = _size[1];
250 std::vector<float> fullresIm(w*h); // accumulate full-resolution image into this buffer
251 const int nbasis = _size.size() > 2 ? _size[2] : 1; // number of basis functions
252
253 /* Sum each component */
254 int const npix = w*h;
255 for (int i = 0; i != nbasis; ++i) {
256 float *pl = &fullresIm[0];
257 float const fac = _poly->basis[i];
258 float const *ppc = &_comp[i*w*h];
259
260 for (int j = 0; j != npix; ++j) {
261 pl[j] += fac*ppc[j];
262 }
263 }
264 /*
265 * We now have the image reconstructed at internal resolution; resample it onto the output resolution
266 * and subpixel offset
267 */
268 float const vigstep = 1/_pixstep;
269 float dx = center[0] - static_cast<int>(center[0]);
270 float dy = center[1] - static_cast<int>(center[1]);
271 if (dx > 0.5) dx -= 1.0;
272 if (dy > 0.5) dy -= 1.0;
273 //
274 // And copy it into place
275 //
276 geom::Box2I bbox = _doComputeBBox(position, center);
277 std::shared_ptr<afw::detection::Psf::Image> im = std::make_shared<afw::detection::Psf::Image>(bbox);
278
279 int sampleW = bbox.getWidth();
280 int sampleH = bbox.getHeight();
281
282 std::vector<float> sampledIm(sampleW*sampleH);
283
284 vignet_resample(&fullresIm[0], w, h,
285 &sampledIm[0], sampleW, sampleH,
286 -dx*vigstep, -dy*vigstep, vigstep, 1.0);
287 {
288 float *pl = &sampledIm[0];
289 float const sum = std::accumulate(pl, pl + sampleW*sampleH, static_cast<float>(0));
290 for (int y = 0; y != sampleH; ++y) {
291 for (int x = 0; x != sampleW; ++x) {
292 (*im)(x, y) = *pl++/sum;
293 }
294 }
295 }
296
297 return im;
298}
299
300/************************************************************************************************************/
301/*
302 * All the rest of this file handles persistence to FITS files
303 */
304namespace table = afw::table;
305
306namespace {
307
308class PsfexPsfSchema1 {
309public:
310 PsfexPsfSchema1() :
311 schema(),
312 ndim(schema.addField<int>("ndim", "Number of elements in group")),
313 ngroup(schema.addField<int>("ngroup", "Number of elements in degree")),
314 ncoeff(schema.addField<int>("ncoeff", "Number of coefficients")),
315
316 _size_size(schema.addField<int>("_size_size", "Size of _size array")),
317 _comp_size(schema.addField<int>("_comp_size", "Size of _comp array")),
318 _context_size(schema.addField<int>("_context_size", "Size of _context array")),
319 // Other scalars
320 averagePosition(afw::table::PointKey<double>::addFields(schema,"averagePosition","average position of stars used to make the PSF","pixel")),
321 _pixstep(schema.addField<float>("_pixstep", "oversampling", "pixel"))
322 {
323 ;
324 }
325
326 table::Schema schema;
327
328 // Sizes in _poly
329 table::Key<int> ndim;
330 table::Key<int> ngroup;
331 table::Key<int> ncoeff;
332 // Sizes of vectors
333 table::Key<int> _size_size;
334 table::Key<int> _comp_size;
335 table::Key<int> _context_size;
336 // Other scalars
337 table::PointKey<double> averagePosition;
338 table::Key<float> _pixstep;
339};
340
341class PsfexPsfSchema2 {
342public:
343 PsfexPsfSchema2(int const ndim, int const ngroup, int const ncoeff,
344 int size_size, int comp_size, int context_size) :
345
346 schema(),
347 group(schema.addField<table::Array<int> >("group", "Groups (of coefficients?)", ndim)),
348 degree(schema.addField<table::Array<int> >("degree", "Degree in each group", ngroup)),
349 basis(schema.addField<table::Array<double> >("basis", "Values of the basis functions", ncoeff)),
350 coeff(schema.addField<table::Array<double> >("coeff", "Polynomial coefficients", ncoeff)),
351 _size(schema.addField<table::Array<int> >("_size", "PSF dimensions", size_size)),
352 _comp(schema.addField<table::Array<float> >("_comp", "Complete pixel data", comp_size)),
353 _context_first( schema.addField<table::Array<double> >("_context_first",
354 "Offset to apply to context data", context_size)),
355 _context_second(schema.addField<table::Array<double> >("_context_second",
356 "Scale to apply to context data", context_size))
357 {
358 ;
359 }
360
361 table::Schema schema;
362 // _poly
363 table::Key<table::Array<int> > group; // len(group) == ndim
364 table::Key<table::Array<int> > degree; // len(degree) == ngroup
365 table::Key<table::Array<double> > basis; // len(basis) == ncoeff
366 table::Key<table::Array<double> > coeff; // len(coeff) == ncoeff
367 // vectors
368 table::Key<table::Array<int> > _size;
369 table::Key<table::Array<float> > _comp;
370 table::Key<table::Array<double> > _context_first;
371 table::Key<table::Array<double> > _context_second;
372};
373
374std::string getPsfexPsfPersistenceName() { return "PsfexPsf"; }
375
376} // anonymous
377
378/************************************************************************************************************/
379
380namespace detail { // PsfexPsfFactory needs to be a friend of PsfexPsf
381class PsfexPsfFactory : public table::io::PersistableFactory {
382public:
383
385read(InputArchive const & archive, CatalogVector const & catalogs) const {
386 LSST_ARCHIVE_ASSERT(catalogs.size() == 2u);
387
389
390 int ndim, ngroup, ncoeff;
391 int size_size, comp_size, context_size;
392 {
393 PsfexPsfSchema1 const & keys = PsfexPsfSchema1();
394 LSST_ARCHIVE_ASSERT(catalogs[0].size() == 1u);
395 LSST_ARCHIVE_ASSERT(catalogs[0].getSchema() == keys.schema);
396 table::BaseRecord const & record = catalogs[0].front();
397
398 // fields in _poly
399 ndim = record.get(keys.ndim);
400 ngroup = record.get(keys.ngroup);
401 ncoeff = record.get(keys.ncoeff);
402 // Other scalars
403 result->_averagePosition = record.get(keys.averagePosition);
404 result->_pixstep = record.get(keys._pixstep);
405 // sizes of vectors
406 size_size = record.get(keys._size_size);
407 comp_size = record.get(keys._comp_size);
408 context_size = record.get(keys._context_size);
409 }
410 // Now we can read the data
411 {
412 PsfexPsfSchema2 const keys(ndim, ngroup, ncoeff,
413 size_size, comp_size, context_size);
414
415 LSST_ARCHIVE_ASSERT(catalogs[1].size() == 1u);
416 LSST_ARCHIVE_ASSERT(catalogs[1].getSchema() == keys.schema);
417 table::BaseRecord const & record = catalogs[1].front();
418
419 // _poly
421 {
422 int const *begin = record.getElement(keys.group);
423 group.assign(begin, begin + ndim);
424
425 for (int i = 0; i != ndim; ++i) {
426 ++group[i]; // poly_init subtracts 1 from each element. Sigh.
427 }
428 }
430 {
431 int const *begin = record.getElement(keys.degree);
432 degree.assign(begin, begin + ngroup);
433 }
434 result->_poly = poly_init(&group[0], group.size(), &degree[0], degree.size());
435 LSST_ARCHIVE_ASSERT(result->_poly->ncoeff == ncoeff);
436
437 {
438 double const *begin = record.getElement(keys.basis);
439 std::copy(begin, begin + ncoeff, result->_poly->basis);
440 }
441 {
442 double const *begin = record.getElement(keys.coeff);
443 std::copy(begin, begin + ncoeff, result->_poly->coeff);
444 }
445 // vectors
446 {
447 int const *begin = record.getElement(keys._size);
448 result->_size.assign(begin, begin + size_size);
449 }
450 {
451 float const *begin = record.getElement(keys._comp);
452 result->_comp.assign(begin, begin + comp_size);
453 }
454 {
455 double const *begin1 = record.getElement(keys._context_first);
456 double const *begin2 = record.getElement(keys._context_second);
457 result->_context.resize(context_size);
458 for (int i = 0; i != context_size; ++i) {
459 result->_context[i].first = begin1[i];
460 result->_context[i].second = begin2[i];
461 }
462 }
463 }
464
465 return result;
466}
467
468explicit PsfexPsfFactory(std::string const & name) : table::io::PersistableFactory(name) {}
469
470};
471}
472
473namespace {
474 detail::PsfexPsfFactory registration(getPsfexPsfPersistenceName());
475}
476
477/************************************************************************************************************/
478
479std::string PsfexPsf::getPythonModule() const { return "lsst.meas.extensions.psfex"; }
480
481std::string PsfexPsf::getPersistenceName() const { return getPsfexPsfPersistenceName(); }
482
484 // First write the dimensions of the various arrays to an HDU, so we can construct the second
485 // HDU using them when we read the Psf back
486 {
487 PsfexPsfSchema1 const keys;
488 afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
490
491 // Sizes in _poly
492 record->set(keys.ndim, _poly->ndim);
493 record->set(keys.ngroup, _poly->ngroup);
494 record->set(keys.ncoeff, _poly->ncoeff);
495 // Other scalars
496 record->set(keys.averagePosition, _averagePosition);
497 record->set(keys._pixstep, _pixstep);
498 record->set(keys._size_size, _size.size());
499 record->set(keys._comp_size, _comp.size());
500 record->set(keys._context_size, _context.size());
501
502 handle.saveCatalog(cat);
503 }
504 // Now we can write the data
505 {
506 PsfexPsfSchema2 const keys(_poly->ndim, _poly->ngroup, _poly->ncoeff,
507 _size.size(), _comp.size(), _context.size());
508 afw::table::BaseCatalog cat = handle.makeCatalog(keys.schema);
509 // _poly
511 {
512 int *begin = record->getElement(keys.group);
513 std::copy(_poly->group, _poly->group + _poly->ndim, begin);
514 }
515 {
516 int *begin = record->getElement(keys.degree);
517 std::copy(_poly->degree, _poly->degree + _poly->ngroup, begin);
518 }
519 {
520 double *begin = record->getElement(keys.basis);
521 std::copy(_poly->basis, _poly->basis + _poly->ncoeff, begin);
522 }
523 {
524 double *begin = record->getElement(keys.coeff);
525 std::copy(_poly->coeff, _poly->coeff + _poly->ncoeff, begin);
526 }
527 // vectors
528 {
529 int *begin = record->getElement(keys._size);
530 std::copy(_size.begin(), _size.end(), begin);
531 }
532 {
533 float *begin = record->getElement(keys._comp);
534 std::copy(_comp.begin(), _comp.end(), begin);
535 }
536 {
537 double *begin1 = record->getElement(keys._context_first);
538 double *begin2 = record->getElement(keys._context_second);
539 for (unsigned int i = 0; i != _context.size(); ++i) {
540 begin1[i] = _context[i].first;
541 begin2[i] = _context[i].second;
542 }
543 }
544
545 handle.saveCatalog(cat);
546 }
547}
548
549}}}}
py::object result
Definition _schema.cc:429
AmpInfoBoxKey bbox
Definition Amplifier.cc:117
#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
table::Schema schema
Definition python.h:134
Basic LSST definitions.
T begin(T... args)
Describe the colour of a source.
Definition Color.h:25
Base class for all records.
Definition BaseRecord.h:31
Field< T >::Value get(Key< T > const &key) const
Return the value of a field for the given key.
Definition BaseRecord.h:151
Field< T >::Element * getElement(Key< T > const &key)
Return a pointer to the underlying elements of a field (non-const).
Definition BaseRecord.h:93
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:489
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.
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:118
void write(lsst::afw::table::io::OutputArchiveHandle &handle) const override
Write the object to one or more catalogs.
Definition PsfexPsf.cc:483
lsst::geom::Box2I doComputeBBox(lsst::geom::Point2D const &position, lsst::afw::image::Color const &color) const override
Compute the bbox of the kernel image at the specified position/color.
Definition PsfexPsf.cc:200
std::shared_ptr< lsst::afw::detection::Psf::Image > doComputeImage(lsst::geom::Point2D const &position, lsst::afw::image::Color const &color) const override
Compute an image of the Psf at the specified position/colour, at pixel position in the output image.
Definition PsfexPsf.cc:188
int getNdim() const
Return the number of dependency parameters in the psfex polynomial fit.
Definition PsfexPsf.cc:113
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:229
std::shared_ptr< lsst::afw::detection::Psf::Image > doComputeKernelImage(lsst::geom::Point2D const &position, lsst::afw::image::Color const &color) const override
Compute an image of the Psf at the specified position/colour, at pixel (0.0, 0.0) in the output image...
Definition PsfexPsf.cc:194
std::string getPersistenceName() const override
Name used in table persistence.
Definition PsfexPsf.cc:481
std::shared_ptr< afw::detection::Psf > resized(int width, int height) const override
Return a clone with specified kernel dimensions.
Definition PsfexPsf.cc:109
std::string getPythonModule() const override
The python module name (for use in table persistence)
Definition PsfexPsf.cc:479
virtual std::shared_ptr< table::io::Persistable > read(InputArchive const &archive, CatalogVector const &catalogs) const
Definition PsfexPsf.cc:385
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++.
T push_back(T... args)
T reserve(T... args)
T size(T... args)
afw::table::PointKey< double > averagePosition
Definition CoaddPsf.cc:354
double w
Definition CoaddPsf.cc:70
table::Key< int > ncoeff
Definition PsfexPsf.cc:331
table::Key< int > _context_size
Definition PsfexPsf.cc:335
table::Key< int > ndim
Definition PsfexPsf.cc:329
table::Key< int > ngroup
Definition PsfexPsf.cc:330
table::Key< table::Array< int > > degree
Definition PsfexPsf.cc:364
table::Key< table::Array< float > > _comp
Definition PsfexPsf.cc:369
table::Key< float > _pixstep
Definition PsfexPsf.cc:338
table::Key< table::Array< double > > _context_first
Definition PsfexPsf.cc:370
table::Key< table::Array< double > > coeff
Definition PsfexPsf.cc:366
table::Key< table::Array< double > > basis
Definition PsfexPsf.cc:365
table::Key< table::Array< double > > _context_second
Definition PsfexPsf.cc:371
table::Key< table::Array< int > > group
Definition PsfexPsf.cc:363
table::Key< int > _comp_size
Definition PsfexPsf.cc:334
table::Key< int > _size_size
Definition PsfexPsf.cc:333
Key< int > psf
Definition Exposure.cc:65