LSST Applications  21.0.0+04719a4bac,21.0.0-1-ga51b5d4+f5e6047307,21.0.0-11-g2b59f77+a9c1acf22d,21.0.0-11-ga42c5b2+86977b0b17,21.0.0-12-gf4ce030+76814010d2,21.0.0-13-g1721dae+760e7a6536,21.0.0-13-g3a573fe+768d78a30a,21.0.0-15-g5a7caf0+f21cbc5713,21.0.0-16-g0fb55c1+b60e2d390c,21.0.0-19-g4cded4ca+71a93a33c0,21.0.0-2-g103fe59+bb20972958,21.0.0-2-g45278ab+04719a4bac,21.0.0-2-g5242d73+3ad5d60fb1,21.0.0-2-g7f82c8f+8babb168e8,21.0.0-2-g8f08a60+06509c8b61,21.0.0-2-g8faa9b5+616205b9df,21.0.0-2-ga326454+8babb168e8,21.0.0-2-gde069b7+5e4aea9c2f,21.0.0-2-gecfae73+1d3a86e577,21.0.0-2-gfc62afb+3ad5d60fb1,21.0.0-25-g1d57be3cd+e73869a214,21.0.0-3-g357aad2+ed88757d29,21.0.0-3-g4a4ce7f+3ad5d60fb1,21.0.0-3-g4be5c26+3ad5d60fb1,21.0.0-3-g65f322c+e0b24896a3,21.0.0-3-g7d9da8d+616205b9df,21.0.0-3-ge02ed75+a9c1acf22d,21.0.0-4-g591bb35+a9c1acf22d,21.0.0-4-g65b4814+b60e2d390c,21.0.0-4-gccdca77+0de219a2bc,21.0.0-4-ge8a399c+6c55c39e83,21.0.0-5-gd00fb1e+05fce91b99,21.0.0-6-gc675373+3ad5d60fb1,21.0.0-64-g1122c245+4fb2b8f86e,21.0.0-7-g04766d7+cd19d05db2,21.0.0-7-gdf92d54+04719a4bac,21.0.0-8-g5674e7b+d1bd76f71f,master-gac4afde19b+a9c1acf22d,w.2021.13
LSST Data Management Base Package
Classes | Public Member Functions | List of all members
lsst::afw::geom::SipApproximation Class Referencefinal

A fitter and results class for approximating a general Transform in a form compatible with FITS WCS persistence. More...

#include <SipApproximation.h>

Classes

struct  Grid
 
struct  Solution
 

Public Member Functions

 SipApproximation (std::shared_ptr< TransformPoint2ToPoint2 > pixelToIwc, lsst::geom::Point2D const &crpix, Eigen::Matrix2d const &cd, lsst::geom::Box2D const &bbox, lsst::geom::Extent2I const &gridShape, int order, bool useInverse=true, double svdThreshold=-1)
 Construct a new approximation by fitting on a grid of points. More...
 
 SipApproximation (std::shared_ptr< TransformPoint2ToPoint2 > pixelToIwc, lsst::geom::Point2D const &crpix, Eigen::Matrix2d const &cd, lsst::geom::Box2D const &bbox, lsst::geom::Extent2I const &gridShape, ndarray::Array< double const, 2 > const &a, ndarray::Array< double const, 2 > const &b, ndarray::Array< double const, 2 > const &ap, ndarray::Array< double const, 2 > const &bp, bool useInverse=true)
 Construct from existing SIP coefficients. More...
 
 SipApproximation (SipApproximation const &)=delete
 
SipApproximationoperator= (SipApproximation const &)=delete
 
 SipApproximation (SipApproximation &&) noexcept=default
 
SipApproximationoperator= (SipApproximation &&) noexcept=default
 
 ~SipApproximation () noexcept
 
int getOrder () const noexcept
 Return the polynomial order of the current solution (same for forward and reverse). More...
 
double getA (int p, int q) const
 Return a coefficient of the forward transform polynomial. More...
 
double getB (int p, int q) const
 Return a coefficient of the forward transform polynomial. More...
 
double getAP (int p, int q) const
 Return a coefficient of the reverse transform polynomial. More...
 
double getBP (int p, int q) const
 Return a coefficient of the reverse transform polynomial. More...
 
Eigen::MatrixXd getA () const noexcept
 Return the coefficients of the forward transform polynomial. More...
 
Eigen::MatrixXd getB () const noexcept
 Return the coefficients of the forward transform polynomial. More...
 
Eigen::MatrixXd getAP () const noexcept
 Return the coefficients of the reverse transform polynomial. More...
 
Eigen::MatrixXd getBP () const noexcept
 Return the coefficients of the reverse transform polynomial. More...
 
lsst::geom::Point2D applyForward (lsst::geom::Point2D const &pix) const
 Convert a point from pixels to intermediate world coordinates. More...
 
std::vector< lsst::geom::Point2DapplyForward (std::vector< lsst::geom::Point2D > const &pix) const
 Convert an array of points from pixels to intermediate world coordinates. More...
 
lsst::geom::Point2D applyInverse (lsst::geom::Point2D const &iwcs) const
 Convert a point from intermediate world coordinates to pixels. More...
 
std::vector< lsst::geom::Point2DapplyInverse (std::vector< lsst::geom::Point2D > const &iwcs) const
 Convert an array of points from intermediate world coordinates to pixels. More...
 
lsst::geom::Extent2D getGridStep () const noexcept
 Return the distance between grid points in pixels. More...
 
lsst::geom::Extent2I getGridShape () const noexcept
 Return the number of grid points in x and y. More...
 
lsst::geom::Box2D getBBox () const noexcept
 Return the pixel-coordinate bounding box over which the approximation should be valid. More...
 
lsst::geom::Point2D getPixelOrigin () const noexcept
 Return the pixel origin of the WCS being approximated. More...
 
Eigen::Matrix2d getCdMatrix () const noexcept
 Return the CD matrix of the WCS being approximated. More...
 
void updateGrid (lsst::geom::Extent2I const &shape)
 Update the grid to the given number of points in x and y. More...
 
void refineGrid (int factor=2)
 Update the grid by making it finer by a given integer factor. More...
 
void fit (int order, double svdThreshold=-1)
 Obtain a new solution at the given order with the current grid. More...
 
std::pair< double, double > computeMaxDeviation () const noexcept
 Return the maximum deviation of the solution from the exact transform on the current grid. More...
 

Detailed Description

A fitter and results class for approximating a general Transform in a form compatible with FITS WCS persistence.

The Simple Imaging Polynomial (SIP) convention (Shupe et al 2005) adds forward and reverse polynomial mappings to a standard projection FITS WCS projection (e.g. "TAN" for gnomonic) that relate Intermediate World Coordinates (see Calabretta & Greisen 2002) to image pixel coordinates. The SIP "forward" transform is defined by polynomial coeffients \(A\) and \(B\) that map pixel coordinates \((u, v)\) to Intermediate World Coordinates \((x, y)\) via

\[ \boldsymbol{S}\left[\begin{array}{c} x \\ y \end{array}\right] \equiv \left[\begin{array}{c} x_s \\ y_s \end{array}\right] = \left[\begin{array}{c} (u - u_0) + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{A}_{p,q} (u - u_0)^p (v - v_0)^q \\ (v - v_0) + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{B}_{p,q} (u - u_0)^p (v - v_0)^q \end{array}\right] \]

The reverse transform has essentially the same form:

\[ \left[\begin{array}{c} u - u_0 \\ v - v_0 \end{array}\right] = \left[\begin{array}{c} x_s + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{AP}_{p,q} x_s^p y_s^q \\ y_s + \displaystyle\sum_{p,q}^{0 \le p + q \le N} \mathrm{BP}_{p,q} x_s^p y_s^q \end{array}\right] \]

In both cases, \((u_0, v_0)\) is the pixel origin (CRPIX in FITS WCS) and \(\boldsymbol{S}\) is the inverse of the Jacobian "CD" matrix. Both CRPIX and CD are considered fixed inputs, and we do not attempt to null the zeroth- and first-order terms of \(A\) and \(B\) (as some SIP fitters do); together, these conventions make solving for the coefficients a much simpler linear problem.

Note
In the implementation, we typically refer to \((u-u_0, v-v_0)\) as dpix (for "pixel delta"), and \((x_s, y_s)\) as siwc (for "scaled intermediate world coordinates").

While LSST WCSs are in general too complex to be described exactly in FITS WCS, they can generally be closely approximated by standard FITS WCS projection with additional SIP distortions. This class fits such an approximation, given a TransformPoint2ToPoint2 object that represents the exact mapping from pixels to Intermediate World Coordinates with a SIP distortion.

Definition at line 94 of file SipApproximation.h.

Constructor & Destructor Documentation

◆ SipApproximation() [1/4]

lsst::afw::geom::SipApproximation::SipApproximation ( std::shared_ptr< TransformPoint2ToPoint2 pixelToIwc,
lsst::geom::Point2D const &  crpix,
Eigen::Matrix2d const &  cd,
lsst::geom::Box2D const &  bbox,
lsst::geom::Extent2I const &  gridShape,
int  order,
bool  useInverse = true,
double  svdThreshold = -1 
)

Construct a new approximation by fitting on a grid of points.

Parameters
[in]pixelToIwcThe true Transform to approximate. Should go from pixels to Intermediate World Coordinates when applyForward is called.
[in]crpixPixel origin, using the LSST 0-indexed convention rather than the FITS 1-indexed convention; equal to (CRPIX1 - 1, CRPIX2 - 1).
[in]cdNominal Jacobian ("CD" in FITS WCS).
[in]bboxPixel-coordinate bounding box over which the approximation should be valid. Used to construct the grid of points to fit.
[in]gridShapeNumber of points in x and y for the grid of points.
[in]orderOrder of the polynomial (same for forward and reverse transforms).
[in]useInverseIf true, the inverse SIP transform will be fit and compared to data points generated by calls to pixelToIwc.applyInverse instead of pixelToIwc.applyForward.
[in]svdThresholdFraction of the largest singular value at which to declare smaller singular values zero in the least squares solution. Negative values use Eigen's internal default.
Exceptions
lsst::pex::exceptions::InvalidParameterErrorThrown if order is negative or gridShape is non-positive.
Exception Safety
strong

Definition at line 210 of file SipApproximation.cc.

219  :
220  _useInverse(useInverse),
221  _pixelToIwc(std::move(pixelToIwc)),
222  _bbox(bbox),
223  _crpix(crpix),
224  _cdInv(lsst::geom::LinearTransform(cd).inverted()),
225  _grid(new Grid(gridShape, *this)),
226  _solution(Solution::fit(order, svdThreshold, *this))
227 {}
AmpInfoBoxKey bbox
Definition: Amplifier.cc:117
table::PointKey< double > crpix
Definition: OldWcs.cc:131
table::Key< table::Array< double > > cd
Definition: OldWcs.cc:132
A 2D linear coordinate transformation.
T move(T... args)
static std::unique_ptr< Solution > fit(int order_, double svdThreshold, SipApproximation const &parent)

◆ SipApproximation() [2/4]

lsst::afw::geom::SipApproximation::SipApproximation ( std::shared_ptr< TransformPoint2ToPoint2 pixelToIwc,
lsst::geom::Point2D const &  crpix,
Eigen::Matrix2d const &  cd,
lsst::geom::Box2D const &  bbox,
lsst::geom::Extent2I const &  gridShape,
ndarray::Array< double const, 2 > const &  a,
ndarray::Array< double const, 2 > const &  b,
ndarray::Array< double const, 2 > const &  ap,
ndarray::Array< double const, 2 > const &  bp,
bool  useInverse = true 
)

Construct from existing SIP coefficients.

This constructor is primarily intended for testing purposes.

Parameters
[in]pixelToIwcThe true Transform to approximate. Should go from pixels to Intermediate World Coordinates when applyForward is called.
[in]crpixPixel origin, using the LSST 0-indexed convention rather than the FITS 1-indexed convention; equal to (CRPIX1 - 1, CRPIX - 1).
[in]cdNominal Jacobian ("CD" in FITS WCS).
[in]bboxPixel-coordinate bounding box over which the approximation should be valid. Used to construct the grid of points to fit.
[in]gridShapeNumber of points in x and y for the grid of points.
[in]aMatrix of A coefficients, with the first dimension corresponding to powers of \((u - u_0)\) and the second corresponding to powers of \((v - v_0)\).
[in]bMatrix of B coefficients, with the first dimension corresponding to powers of \((u - u_0)\) and the second corresponding to powers of \((v - v_0)\).
[in]apMatrix of AP coefficients, with the first dimension corresponding to powers of \(x_s\) and the second corresponding to powers of \(y_s\).
[in]bpMatrix of BP coefficients, with the first dimension corresponding to powers of \(x_s\) and the second corresponding to powers of \(y_s\).
[in]useInverseIf true, the inverse SIP transform will be compared to data points generated by calls to pixelToIwc.applyInverse instead of pixelToIwc.applyForward.
Exceptions
lsst::pex::exceptions::InvalidParameterErrorThrown if gridShape is non-positive, or any matrix argument is non-square.
Exception Safety
strong

Definition at line 229 of file SipApproximation.cc.

240  :
241  _useInverse(useInverse),
242  _pixelToIwc(std::move(pixelToIwc)),
243  _bbox(bbox),
244  _crpix(crpix),
245  _cdInv(lsst::geom::LinearTransform(cd).inverted()),
246  _grid(new Grid(gridShape, *this)),
247  _solution(
248  new Solution(
249  makePolynomialFromCoeffMatrix(a),
250  makePolynomialFromCoeffMatrix(b),
251  makePolynomialFromCoeffMatrix(ap),
252  makePolynomialFromCoeffMatrix(bp)
253  )
254  )
255 {}
table::Key< int > b
table::Key< int > a

◆ SipApproximation() [3/4]

lsst::afw::geom::SipApproximation::SipApproximation ( SipApproximation const &  )
delete

◆ SipApproximation() [4/4]

lsst::afw::geom::SipApproximation::SipApproximation ( SipApproximation &&  )
defaultnoexcept

◆ ~SipApproximation()

lsst::afw::geom::SipApproximation::~SipApproximation ( )
noexcept

Definition at line 257 of file SipApproximation.cc.

257 {}

Member Function Documentation

◆ applyForward() [1/2]

lsst::geom::Point2D lsst::afw::geom::SipApproximation::applyForward ( lsst::geom::Point2D const &  pix) const

Convert a point from pixels to intermediate world coordinates.

This method is inefficient and should only be used for diagnostic purposes.

Exception Safety
strong

Definition at line 322 of file SipApproximation.cc.

322  {
323  auto cd = _cdInv.inverted();
324  auto ws = _solution->makeWorkspace();
325  return cd(_solution->applyForward(pix - _crpix, ws));
326 }
LinearTransform const inverted() const
Return the inverse transform.
T ws(T... args)

◆ applyForward() [2/2]

std::vector< lsst::geom::Point2D > lsst::afw::geom::SipApproximation::applyForward ( std::vector< lsst::geom::Point2D > const &  pix) const

Convert an array of points from pixels to intermediate world coordinates.

Exception Safety
strong

Definition at line 328 of file SipApproximation.cc.

329  {
330  auto ws = _solution->makeWorkspace();
332  iwc.reserve(pix.size());
333  auto cd = _cdInv.inverted();
334  for (auto const & point : pix) {
335  iwc.push_back(cd(_solution->applyForward(point - _crpix, ws)));
336  }
337  return iwc;
338 }
T push_back(T... args)
T reserve(T... args)
T size(T... args)

◆ applyInverse() [1/2]

lsst::geom::Point2D lsst::afw::geom::SipApproximation::applyInverse ( lsst::geom::Point2D const &  iwcs) const

Convert a point from intermediate world coordinates to pixels.

This method is inefficient and should only be used for diagnostic purposes.

Exception Safety
strong

Definition at line 340 of file SipApproximation.cc.

340  {
341  auto ws = _solution->makeWorkspace();
342  return _solution->applyInverse(_cdInv(iwc), ws) + _crpix;
343 }

◆ applyInverse() [2/2]

std::vector< lsst::geom::Point2D > lsst::afw::geom::SipApproximation::applyInverse ( std::vector< lsst::geom::Point2D > const &  iwcs) const

Convert an array of points from intermediate world coordinates to pixels.

Exception Safety
strong

Definition at line 345 of file SipApproximation.cc.

346  {
347  auto ws = _solution->makeWorkspace();
349  pix.reserve(iwc.size());
350  for (auto const & point : iwc) {
351  pix.push_back(_solution->applyInverse(_cdInv(point), ws) + _crpix);
352  }
353  return pix;
354 }

◆ computeMaxDeviation()

std::pair< double, double > lsst::afw::geom::SipApproximation::computeMaxDeviation ( ) const
noexcept

Return the maximum deviation of the solution from the exact transform on the current grid.

The deviations are in scaled intermediate world coordinates \(\sqrt{\delta x_s^2 \delta y_s^2}\) for the forward transform and in pixels \((\delta u^2, \delta v^2)\) for the reverse transform (respectively). Note that in the common case where the CD matrix includes the scaling from angle units to pixel units, the scaled intermediate world coordinate values are also in (nominal) pixel units.

Definition at line 381 of file SipApproximation.cc.

381  {
382  std::pair<double, double> maxDiff(0.0, 0.0);
383  auto ws = _solution->makeWorkspace();
384  for (std::size_t i = 0; i < _grid->dpix1.size(); ++i) {
385  auto siwc2 = _solution->applyForward(_grid->dpix1[i], ws);
386  auto dpix2 = _solution->applyInverse(_grid->siwc[i], ws);
387  maxDiff.first = std::max(maxDiff.first, (_grid->siwc[i] - siwc2).computeNorm());
388  maxDiff.second = std::max(maxDiff.second, (_grid->dpix2[i] - dpix2).computeNorm());
389  }
390  return maxDiff;
391 }
T max(T... args)

◆ fit()

void lsst::afw::geom::SipApproximation::fit ( int  order,
double  svdThreshold = -1 
)

Obtain a new solution at the given order with the current grid.

Parameters
[in]orderPolynomial order to fit.
[in]svdThresholdFraction of the largest singular value at which to declare smaller singular values zero in the least squares solution. Negative values use Eigen's internal default.
Exceptions
pex::exceptions::LogicErrorThrown if the number of free parameters implied by order is larger than the number of data points defined by the grid.
Exception Safety
strong

Definition at line 377 of file SipApproximation.cc.

377  {
378  _solution = Solution::fit(order, svdThreshold, *this);
379 }

◆ getA() [1/2]

Eigen::MatrixXd lsst::afw::geom::SipApproximation::getA ( ) const
noexcept

Return the coefficients of the forward transform polynomial.

Definition at line 294 of file SipApproximation.cc.

294  {
295  return makeCoefficientMatrix(
296  getOrder(),
297  [this](int p, int q) { return getA(p, q); }
298  );
299 }
Eigen::MatrixXd getA() const noexcept
Return the coefficients of the forward transform polynomial.
int getOrder() const noexcept
Return the polynomial order of the current solution (same for forward and reverse).

◆ getA() [2/2]

double lsst::afw::geom::SipApproximation::getA ( int  p,
int  q 
) const

Return a coefficient of the forward transform polynomial.

Out-of-bounds arguments yields undefined behavior.

Exception Safety
strong

Definition at line 263 of file SipApproximation.cc.

263  {
264  return _solution->a[_solution->a.getBasis().index(p, q)];
265 }

◆ getAP() [1/2]

Eigen::MatrixXd lsst::afw::geom::SipApproximation::getAP ( ) const
noexcept

Return the coefficients of the reverse transform polynomial.

Definition at line 308 of file SipApproximation.cc.

308  {
309  return makeCoefficientMatrix(
310  getOrder(),
311  [this](int p, int q) { return getAP(p, q); }
312  );
313 }
Eigen::MatrixXd getAP() const noexcept
Return the coefficients of the reverse transform polynomial.

◆ getAP() [2/2]

double lsst::afw::geom::SipApproximation::getAP ( int  p,
int  q 
) const

Return a coefficient of the reverse transform polynomial.

Out-of-bounds arguments yields undefined behavior.

Exception Safety
strong

Definition at line 271 of file SipApproximation.cc.

271  {
272  return _solution->ap[_solution->ap.getBasis().index(p, q)];
273 }

◆ getB() [1/2]

Eigen::MatrixXd lsst::afw::geom::SipApproximation::getB ( ) const
noexcept

Return the coefficients of the forward transform polynomial.

Definition at line 301 of file SipApproximation.cc.

301  {
302  return makeCoefficientMatrix(
303  getOrder(),
304  [this](int p, int q) { return getB(p, q); }
305  );
306 }
Eigen::MatrixXd getB() const noexcept
Return the coefficients of the forward transform polynomial.

◆ getB() [2/2]

double lsst::afw::geom::SipApproximation::getB ( int  p,
int  q 
) const

Return a coefficient of the forward transform polynomial.

Out-of-bounds arguments yields undefined behavior.

Exception Safety
strong

Definition at line 267 of file SipApproximation.cc.

267  {
268  return _solution->b[_solution->b.getBasis().index(p, q)];
269 }

◆ getBBox()

lsst::geom::Box2D lsst::afw::geom::SipApproximation::getBBox ( ) const
inlinenoexcept

Return the pixel-coordinate bounding box over which the approximation should be valid.

Definition at line 291 of file SipApproximation.h.

291 { return _bbox; }

◆ getBP() [1/2]

Eigen::MatrixXd lsst::afw::geom::SipApproximation::getBP ( ) const
noexcept

Return the coefficients of the reverse transform polynomial.

Definition at line 315 of file SipApproximation.cc.

315  {
316  return makeCoefficientMatrix(
317  getOrder(),
318  [this](int p, int q) { return getBP(p, q); }
319  );
320 }
Eigen::MatrixXd getBP() const noexcept
Return the coefficients of the reverse transform polynomial.

◆ getBP() [2/2]

double lsst::afw::geom::SipApproximation::getBP ( int  p,
int  q 
) const

Return a coefficient of the reverse transform polynomial.

Out-of-bounds arguments yields undefined behavior.

Exception Safety
strong

Definition at line 275 of file SipApproximation.cc.

275  {
276  return _solution->bp[_solution->bp.getBasis().index(p, q)];
277 }

◆ getCdMatrix()

Eigen::Matrix2d lsst::afw::geom::SipApproximation::getCdMatrix ( ) const
inlinenoexcept

Return the CD matrix of the WCS being approximated.

Definition at line 297 of file SipApproximation.h.

297 { return _cdInv.inverted().getMatrix(); }
Matrix const & getMatrix() const noexcept

◆ getGridShape()

lsst::geom::Extent2I lsst::afw::geom::SipApproximation::getGridShape ( ) const
noexcept

Return the number of grid points in x and y.

Definition at line 361 of file SipApproximation.cc.

361  {
362  return _grid->shape;
363 }

◆ getGridStep()

lsst::geom::Extent2D lsst::afw::geom::SipApproximation::getGridStep ( ) const
noexcept

Return the distance between grid points in pixels.

Definition at line 356 of file SipApproximation.cc.

356  {
357  return lsst::geom::Extent2D(_bbox.getWidth()/_grid->shape.getX(),
358  _bbox.getHeight()/_grid->shape.getY());
359 }
double getWidth() const noexcept
1-d interval accessors
Definition: Box.h:529
double getHeight() const noexcept
1-d interval accessors
Definition: Box.h:530
Extent< double, 2 > Extent2D
Definition: Extent.h:400

◆ getOrder()

int lsst::afw::geom::SipApproximation::getOrder ( ) const
noexcept

Return the polynomial order of the current solution (same for forward and reverse).

Definition at line 259 of file SipApproximation.cc.

259  {
260  return _solution->a.getBasis().getOrder();
261 }

◆ getPixelOrigin()

lsst::geom::Point2D lsst::afw::geom::SipApproximation::getPixelOrigin ( ) const
inlinenoexcept

Return the pixel origin of the WCS being approximated.

Definition at line 294 of file SipApproximation.h.

294 { return lsst::geom::Point2D(_crpix); }
Point< double, 2 > Point2D
Definition: Point.h:324

◆ operator=() [1/2]

SipApproximation& lsst::afw::geom::SipApproximation::operator= ( SipApproximation &&  )
defaultnoexcept

◆ operator=() [2/2]

SipApproximation& lsst::afw::geom::SipApproximation::operator= ( SipApproximation const &  )
delete

◆ refineGrid()

void lsst::afw::geom::SipApproximation::refineGrid ( int  factor = 2)

Update the grid by making it finer by a given integer factor.

Exceptions
lsst::pex::exceptions::InvalidParameterErrorThrown if factor is non-positive.
Exception Safety
strong

Definition at line 369 of file SipApproximation.cc.

369  {
370  // We shrink the grid spacing by the given factor, which is not the same
371  // as increasing the number of grid points by that factor, because there
372  // is one more grid point that step in each dimension.
373  lsst::geom::Extent2I unit(1);
374  updateGrid((_grid->shape - unit)*f + unit);
375 }
void updateGrid(lsst::geom::Extent2I const &shape)
Update the grid to the given number of points in x and y.

◆ updateGrid()

void lsst::afw::geom::SipApproximation::updateGrid ( lsst::geom::Extent2I const &  shape)

Update the grid to the given number of points in x and y.

This does not invalidate or modify the current solution; this allows the user to fit with a coarse grid and then check whether the solution still works well on a finer grid.

Exceptions
lsst::pex::exceptions::InvalidParameterErrorThrown if shape is non-positive.
Exception Safety
strong

Definition at line 365 of file SipApproximation.cc.

365  {
366  _grid = std::make_unique<Grid>(shape, *this);
367 }

The documentation for this class was generated from the following files: