File AstrometryTransform.h¶
- 
namespace 
lsst Class for a simple mapping implementing a generic AstrometryTransform.
Remove all non-astronomical counts from the Chunk Exposure’s pixels.
Forward declarations for lsst::utils::Cache
For details on the Cache class, see the Cache.h file.
It uses a template rather than a pointer so that the derived classes can use the specifics of the transform. The class simplePolyMapping overloads a few routines.
A base class for image defects
Numeric constants used by the Integrate.h integrator routines.
Compute Image Statistics
- Note
 Gauss-Kronrod-Patterson quadrature coefficients for use in quadpack routine qng. These coefficients were calculated with 101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov 1981.
- Note
 The Statistics class itself can only handle lsst::afw::image::MaskedImage() types. The philosophy has been to handle other types by making them look like lsst::afw::image::MaskedImage() and reusing that code. Users should have no need to instantiate a Statistics object directly, but should use the overloaded makeStatistics() factory functions.
- 
namespace 
jointcal Typedefs
- 
typedef 
void() lsst::jointcal::AstrometryTransformFun(const double, const double, double &, double &, const void *) signature of the user-provided routine that actually does the coordinate transform for UserTransform.
Functions
- 
std::ostream &
operator<<(std::ostream &stream, AstrometryTransform const &transform)¶ Delegates to transform.dump()
- 
std::unique_ptr<AstrometryTransform> 
compose(AstrometryTransform const &left, AstrometryTransform const &right)¶ Returns a pointer to a composition of transforms, representing
left(right()).Deletion of returned value to be done by caller.
If
left->composeAndReduce(right)returns NULL, build a AstrometryTransformComposition and return it. This routine implements “run-time” compositions. When there is a possible “reduction” (e.g. compositions of polynomials), compose detects it and returns a genuine AstrometryTransform.- Return
 The composed transform.
- 
std::unique_ptr<AstrometryTransform> 
compose(AstrometryTransform const &left, AstrometryTransformIdentity const &right)¶ 
- 
bool 
isIntegerShift(const AstrometryTransform *transform)¶ Shorthand test to tell if a transform is a simple integer shift.
- 
std::shared_ptr<AstrometryTransformPolynomial> 
inversePolyTransform(AstrometryTransform const &forward, Frame const &domain, double const precision, std::size_t maxOrder = 9, std::size_t nSteps = 50)¶ Approximate the inverse by a polynomial, to some precision.
- Return
 A polynomial that best approximates forward.
- Parameters
 forward: Transform to be inverted.[in] domain: The domain of forward.[in] precision: Require that \(chi2 / (nsteps^2) < precision^2\).[in] maxOrder: The maximum order allowed of the inverse polynomial.[in] nSteps: The number of sample points per axis (nSteps^2 total points).
- 
AstrometryTransformLinear 
normalizeCoordinatesTransform(const Frame &frame)¶ 
- 
std::unique_ptr<AstrometryTransform> 
astrometryTransformRead(const std::string &fileName)¶ The virtual constructor from a file.
- 
std::unique_ptr<AstrometryTransform> 
astrometryTransformRead(std::istream &s)¶ The virtual constructor from a file.
- 
class 
AstrometryTransform - #include <AstrometryTransform.h>
a virtual (interface) class for geometric transformations.
We implement here One AstrometryTransform interface class, and actual derived classes. Composition in the usual (mathematical) sense is provided using compose(), and some classes (e.g. AstrometryTransformLinear
) handle a * operator. Generic inversion by iteration exists, but it is at least 10 times slower than the corresponding “direct
transformation”. If a transform has an analytical inverse, then providing inverseTransform is obviously a very good idea. Before resorting to inverseTransform, consider using
StarMatchList::inverseTransform(). AstrometryTransformLinear::inverted() and TanPixelToRaDec::inverted() exist. The classes also provide derivation and linear approximation.Subclassed by lsst::jointcal::AstrometryTransformIdentity, lsst::jointcal::AstrometryTransformPolynomial, lsst::jointcal::AstrometryTransformSkyWcs, lsst::jointcal::BaseTanWcs, lsst::jointcal::TanRaDecToPixel, lsst::jointcal::UserTransform
Public Functions
- 
virtual void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const = 0 
- 
void 
apply(Point const &in, Point &out) const applies the tranfo to in and writes into out. Is indeed virtual.
- 
Point 
apply(Point const &in) const All these apply(..) shadow the virtual one in derived classes, unless one writes “using
AstrometryTransform::apply”.
- 
Frame 
apply(Frame const &inputframe, bool inscribed) const Transform a bounding box, taking either the inscribed or circumscribed box.
- Return
 The transformed frame.
- Parameters
 [in] inputframe: The frame to be transformed.[in] inscribed: Return the inscribed (true) or circumscribed (false) box.
- 
virtual void 
dump(std::ostream &stream = std::cout) const = 0 dumps the transform coefficients to stream.
- 
std::string 
__str__() 
- 
virtual double 
fit(StarMatchList const &starMatchList) = 0 fits a transform to a std::list of Point pairs (p1,p2, the Point fields in StarMatch).
After the fit this(p1) yields approximately p2. The returned value is the sum of squared residuals. If you want to fit a partial transform (e.g. such that this(T1(p1)) = T2(p2), use StarMatchList::applyTransform beforehand.
- 
void 
transformStar(FatPoint &in) const 
- 
virtual double 
getJacobian(Point const &point) const returns the local jacobian.
- 
virtual std::unique_ptr<AstrometryTransform> 
clone() const = 0 returns a copy (allocated by new) of the transformation.
- 
virtual std::unique_ptr<AstrometryTransform> 
composeAndReduce(AstrometryTransform const &right) const Return a reduced composition of newTransform = this(right()), or nullptr if it cannot be reduced.
“Reduced” in this context means that they are capable of being merged into a single transform, for example, for two polynomials:
\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]we would haveh = f.composeAndReduce(g) == 2 - 6x + 9x^2.To be overloaded by derived classes if they can properly reduce the composition.
- Return
 The new reduced and composed AstrometryTransform, or nullptr if no such reduction is possible.
- Parameters
 right: The transform to apply first.
- 
virtual double 
getJacobian(const double x, const double y) const returns the local jacobian.
- 
virtual void 
computeDerivative(Point const &where, AstrometryTransformLinear &derivative, const double step = 0.01) const Computes the local Derivative of a transform, w.r.t. position.
Step is used for numerical derivation.
- 
virtual AstrometryTransformLinear 
linearApproximation(Point const &where, const double step = 0.01) const linear (local) approximation.
- 
virtual void 
transformErrors(Point const &where, const double *vIn, double *vOut) const transform errors (represented as double[3] in order V(xx),V(yy),Cov(xy))
- 
virtual std::unique_ptr<AstrometryTransform> 
inverseTransform(const double precision, const Frame ®ion) const returns an inverse transform. Numerical if not overloaded.
precision and region refer to the “input” side of this, and hence to the output side of the returned AstrometryTransform.
- 
void 
getParams(double *params) const params should be at least Npar() long
- 
void 
offsetParams(Eigen::VectorXd const &delta) 
- 
virtual double 
paramRef(Eigen::Index const i) const 
- 
virtual double &
paramRef(Eigen::Index const i) 
- 
virtual void 
paramDerivatives(Point const &where, double *dx, double *dy) const Derivative w.r.t parameters. Derivatives should be al least 2*NPar long. first Npar, for x, last Npar for y.
- 
virtual std::unique_ptr<AstrometryTransform> 
roughInverse(const Frame ®ion) const Rough inverse.
Stored by the numerical inverter to guess starting point for the trials. Just here to enable overloading.
- 
virtual std::size_t 
getNpar() const returns the number of parameters (to compute chi2’s)
- 
virtual std::shared_ptr<ast::Mapping> 
toAstMap(jointcal::Frame const &domain) const Create an equivalent AST mapping for this transformation, including an analytic inverse if possible.
- Return
 An AST Mapping that represents this transformation.
- Parameters
 domain: The domain of the transform, to help find an inverse.
- 
void 
write(const std::string &fileName) const 
- 
virtual void 
write(std::ostream &stream) const 
- 
virtual 
~AstrometryTransform() 
 - 
virtual void 
 
- 
class 
AstrometryTransformIdentity: public lsst::jointcal::AstrometryTransform - #include <AstrometryTransform.h>
A do-nothing transformation. It anyway has dummy routines to mimick a AstrometryTransform.
Public Functions
- 
AstrometryTransformIdentity() constructor.
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const xOut = xIn; yOut = yIn !
- 
double 
fit(StarMatchList const &starMatchList) fits a transform to a std::list of Point pairs (p1,p2, the Point fields in StarMatch).
After the fit this(p1) yields approximately p2. The returned value is the sum of squared residuals. If you want to fit a partial transform (e.g. such that this(T1(p1)) = T2(p2), use StarMatchList::applyTransform beforehand.
- 
std::unique_ptr<AstrometryTransform> 
composeAndReduce(AstrometryTransform const &right) const Return a reduced composition of newTransform = this(right()), or nullptr if it cannot be reduced.
“Reduced” in this context means that they are capable of being merged into a single transform, for example, for two polynomials:
\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]we would haveh = f.composeAndReduce(g) == 2 - 6x + 9x^2.To be overloaded by derived classes if they can properly reduce the composition.
- Return
 The new reduced and composed AstrometryTransform, or nullptr if no such reduction is possible.
- Parameters
 right: The transform to apply first.
- 
std::size_t 
getNpar() const returns the number of parameters (to compute chi2’s)
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
void 
computeDerivative(Point const &where, AstrometryTransformLinear &derivative, const double step = 0.01) const Computes the local Derivative of a transform, w.r.t. position.
Step is used for numerical derivation.
- 
virtual AstrometryTransformLinear 
linearApproximation(Point const &where, const double step = 0.01) const linear approximation.
- 
std::shared_ptr<ast::Mapping> 
toAstMap(jointcal::Frame const &domain) const Create an equivalent AST mapping for this transformation, including an analytic inverse if possible.
- Return
 An AST Mapping that represents this transformation.
- Parameters
 domain: The domain of the transform, to help find an inverse.
- 
void 
write(std::ostream &s) const 
- 
void 
read(std::istream &s) 
 - 
 
- 
class 
AstrometryTransformLinear: public lsst::jointcal::AstrometryTransformPolynomial - #include <AstrometryTransform.h>
implements the linear transformations (6 real coefficients).
Subclassed by lsst::jointcal::AstrometryTransformLinearRot, lsst::jointcal::AstrometryTransformLinearScale, lsst::jointcal::AstrometryTransformLinearShift
Public Functions
- 
AstrometryTransformLinear() the default constructor constructs the do-nothing transformation.
- 
AstrometryTransformLinear(AstrometryTransformPolynomial const &transform) This triggers an exception if P.getOrder() != 1.
- 
AstrometryTransformLinear 
operator*(AstrometryTransformLinear const &right) const enables to combine linear tranformations: T1=T2*T3 is legal.
- 
AstrometryTransformLinear 
inverted() const returns the inverse: T1 = T2.inverted();
- 
void 
computeDerivative(Point const &where, AstrometryTransformLinear &derivative, const double step = 0.01) const specialised analytic routine
- 
AstrometryTransformLinear 
linearApproximation(Point const &where, const double step = 0.01) const linear (local) approximation.
- 
AstrometryTransformLinear(const double ox, const double oy, const double aa11, const double aa12, const double aa21, const double aa22) Construct a AstrometryTransformLinear from parameters.
- 
AstrometryTransformLinear(AstrometryTransformIdentity const&) Handy converter:
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
std::unique_ptr<AstrometryTransform> 
inverseTransform(const double precision, const Frame ®ion) const returns an inverse transform. Numerical if not overloaded.
precision and region refer to the “input” side of this, and hence to the output side of the returned AstrometryTransform.
- 
double 
A11() const 
- 
double 
A12() const 
- 
double 
A21() const 
- 
double 
A22() const 
- 
double 
Dx() const 
- 
double 
Dy() const 
Protected Functions
- 
double &
a11() 
- 
double &
a12() 
- 
double &
a21() 
- 
double &
a22() 
- 
double &
dx() 
- 
double &
dy() 
Friends
- 
friend 
lsst::jointcal::AstrometryTransform 
- 
friend 
lsst::jointcal::AstrometryTransformIdentity 
- 
friend 
lsst::jointcal::AstrometryTransformPolynomial 
 - 
 
- 
class 
AstrometryTransformLinearRot: public lsst::jointcal::AstrometryTransformLinear - #include <AstrometryTransform.h>
just here to provide a specialized constructor, and fit.
Public Functions
- 
AstrometryTransformLinearRot() 
- 
AstrometryTransformLinearRot(const double angleRad, const Point *center = nullptr, const double scaleFactor = 1.0) 
- 
double 
fit(StarMatchList const &starMatchList) guess what
- 
std::size_t 
getNpar() const total number of parameters
 - 
 
- 
class 
AstrometryTransformLinearScale: public lsst::jointcal::AstrometryTransformLinear - #include <AstrometryTransform.h>
just here to provide specialized constructors. AstrometryTransformLinear fit routine.
Public Functions
- 
AstrometryTransformLinearScale(const double scale = 1) 
- 
AstrometryTransformLinearScale(const double scaleX, const double scaleY) 
- 
std::size_t 
getNpar() const total number of parameters
 - 
 
- 
class 
AstrometryTransformLinearShift: public lsst::jointcal::AstrometryTransformLinear - #include <AstrometryTransform.h>
just here to provide a specialized constructor, and fit.
Public Functions
- 
AstrometryTransformLinearShift(double ox = 0., double oy = 0.) Add ox and oy.
- 
AstrometryTransformLinearShift(Point const &point) 
- 
double 
fit(StarMatchList const &starMatchList) guess what
- 
std::size_t 
getNpar() const total number of parameters
 - 
 
- 
class 
AstrometryTransformPolynomial: public lsst::jointcal::AstrometryTransform - #include <AstrometryTransform.h>
Polynomial transformation class.
Subclassed by lsst::jointcal::AstrometryTransformLinear
Public Functions
- 
AstrometryTransformPolynomial(std::size_t order = 1) Default transform : identity for all orders (>=1 ).
- Parameters
 order: The highest total power (x+y) of monomials of this polynomial.
- 
AstrometryTransformPolynomial(const AstrometryTransform *transform, const Frame &frame, std::size_t order, std::size_t nPoint = 1000) Constructs a “polynomial image” from an existing transform, over a specified domain.
Constructs a polynomial approximation to an afw::geom::TransformPoint2ToPoint2.
- Parameters
 [in] transform: The transform to be approximated.[in] domain: The valid domain of the transform.[in] order: The polynomial order to use when approximating.[in] nSteps: The number of sample points per axis (nSteps^2 total points).
- 
void 
setOrder(std::size_t order) Sets the polynomial order (the highest sum of exponents of the largest monomial).
- 
std::size_t 
getOrder() const Returns the polynomial order.
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const 
- 
void 
computeDerivative(Point const &where, AstrometryTransformLinear &derivative, const double step = 0.01) const specialised analytic routine
- 
virtual void 
transformPosAndErrors(const FatPoint &in, FatPoint &out) const a mix of apply and Derivative
- 
std::size_t 
getNpar() const total number of parameters
- 
double 
fit(StarMatchList const &starMatchList) guess what
- 
AstrometryTransformPolynomial 
operator*(AstrometryTransformPolynomial const &right) const Composition (internal stuff in quadruple precision)
- 
AstrometryTransformPolynomial 
operator+(AstrometryTransformPolynomial const &right) const Addition.
- 
AstrometryTransformPolynomial 
operator-(AstrometryTransformPolynomial const &right) const Subtraction.
- 
std::unique_ptr<AstrometryTransform> 
composeAndReduce(AstrometryTransformPolynomial const &right) const Return a reduced composition of newTransform = this(right()), or nullptr if it cannot be reduced.
“Reduced” in this context means that they are capable of being merged into a single transform, for example, for two polynomials:
\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]we would haveh = f.composeAndReduce(g) == 2 - 6x + 9x^2.To be overloaded by derived classes if they can properly reduce the composition.
- Return
 The new reduced and composed AstrometryTransform, or nullptr if no such reduction is possible.
- Parameters
 right: The transform to apply first.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
double 
coeff(std::size_t powX, std::size_t powY, std::size_t whichCoord) const access to coefficients (read only)
- 
double 
coeffOrZero(std::size_t powX, std::size_t powY, std::size_t whichCoord) const read access, zero if beyond order
- 
double 
determinant() const 
- 
double 
paramRef(Eigen::Index const i) const 
- 
double &
paramRef(Eigen::Index const i) 
- 
void 
paramDerivatives(Point const &where, double *dx, double *dy) const Derivative w.r.t parameters. Derivatives should be al least 2*NPar long. first Npar, for x, last Npar for y.
- 
std::shared_ptr<ast::Mapping> 
toAstMap(jointcal::Frame const &domain) const Create an equivalent AST mapping for this transformation, including an analytic inverse if possible.
- Return
 An AST Mapping that represents this transformation.
- Parameters
 domain: The domain of the transform, to help find an inverse.
- 
void 
write(std::ostream &s) const 
- 
void 
read(std::istream &s) 
Private Functions
- 
double 
computeFit(StarMatchList const &starMatchList, AstrometryTransform const &shiftToCenter, const bool useErrors)¶ 
- 
void 
computeMonomials(double xIn, double yIn, double *monomial) const¶ 
- 
ndarray::Array<double, 2, 2> 
toAstPolyMapCoefficients() const¶ Return the sparse coefficients matrix that ast::PolyMap requires.
- See
 ast::PolyMap for details of the structure of this matrix.
 - 
 
- 
class 
AstrometryTransformSkyWcs: public lsst::jointcal::AstrometryTransform - #include <AstrometryTransform.h>
A AstrometryTransform that holds a SkyWcs
This is intended to hold the initial estimate for the WCS. It need not be TAN-SIP, nor exactly representable as a FITS WCS.
AstrometryTransformSkyWcs does not inherit from BaseTanWcs for two reasons:
There is no need.
It is not clear how to implement all of the BaseTanWcs interface, especially corrections.
Public Functions
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const 
- 
double 
fit(const StarMatchList &starMatchList) Not implemented; throws pex::exceptions::LogicError.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
 
- 
class 
BaseTanWcs: public lsst::jointcal::AstrometryTransform Subclassed by lsst::jointcal::TanPixelToRaDec, lsst::jointcal::TanSipPixelToRaDec
Public Functions
- 
BaseTanWcs(AstrometryTransformLinear const &pixToTan, Point const &tangentPoint, const AstrometryTransformPolynomial *corrections = nullptr) 
- 
BaseTanWcs(const BaseTanWcs &original) 
- 
void 
operator=(const BaseTanWcs &original) 
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const Transform pixels to ICRS RA, Dec in degrees.
- 
Point 
getTangentPoint() const Get the sky origin (CRVAL in FITS WCS terminology) in degrees.
- 
AstrometryTransformLinear 
getLinPart() const The Linear part (corresponding to CD’s and CRPIX’s)
- 
const AstrometryTransformPolynomial *
getCorr() const Get a non-owning pointer to the correction transform polynomial.
- 
void 
setCorrections(std::unique_ptr<AstrometryTransformPolynomial> corrections) Assign the correction polynomial (what it means is left to derived classes)
- 
Point 
getCrPix() const Get the pixel origin of the WCS (CRPIX in FITS WCS terminology, but zero-based)
- 
virtual AstrometryTransformPolynomial 
getPixelToTangentPlane() const = 0 Get a transform from pixels to tangent plane (degrees) This is a linear transform plus the effects of the correction
- 
virtual void 
pixToTangentPlane(double xPixel, double yPixel, double &xTangentPlane, double &yTangentPlane) const = 0 Transform from pixels to tangent plane (degrees)
- 
~BaseTanWcs() 
Protected Attributes
- 
AstrometryTransformLinear 
linPixelToTan 
- 
std::unique_ptr<AstrometryTransformPolynomial> 
corr 
- 
double 
ra0 
- 
double 
dec0 
- 
double 
cos0 
- 
double 
sin0 
- 
 
- 
class 
TanPixelToRaDec: public lsst::jointcal::BaseTanWcs - #include <AstrometryTransform.h>
The transformation that handles pixels to sideral transformations (Gnomonic, possibly with polynomial distortions).
Public Functions
- 
TanPixelToRaDec(AstrometryTransformLinear const &pixToTan, Point const &tangentPoint, const AstrometryTransformPolynomial *corrections = nullptr) pixToTan describes the transform from pix to tangent plane (degrees). TangentPoint in degrees. Corrections are applied between Lin and deprojection parts (as in Swarp).
- 
AstrometryTransformPolynomial 
getPixelToTangentPlane() const the transformation from pixels to tangent plane (degrees)
- 
virtual void 
pixToTangentPlane(double xPixel, double yPixel, double &xTangentPlane, double &yTangentPlane) const transforms from pixel space to tangent plane (degrees)
- 
TanPixelToRaDec() 
- 
TanPixelToRaDec 
operator*(AstrometryTransformLinear const &right) const composition with AstrometryTransformLinear
- 
std::unique_ptr<AstrometryTransform> 
composeAndReduce(AstrometryTransformLinear const &right) const Return a reduced composition of newTransform = this(right()), or nullptr if it cannot be reduced.
“Reduced” in this context means that they are capable of being merged into a single transform, for example, for two polynomials:
\[ f(x) = 1 + x^2, g(x) = -1 + 3x \]we would haveh = f.composeAndReduce(g) == 2 - 6x + 9x^2.To be overloaded by derived classes if they can properly reduce the composition.
- Return
 The new reduced and composed AstrometryTransform, or nullptr if no such reduction is possible.
- Parameters
 right: The transform to apply first.
- 
TanRaDecToPixel 
inverted() const approximate inverse : it ignores corrections;
- 
std::unique_ptr<AstrometryTransform> 
roughInverse(const Frame ®ion) const Overload the “generic routine” (available for all AstrometryTransform types.
- 
std::unique_ptr<AstrometryTransform> 
inverseTransform(const double precision, const Frame ®ion) const Inverse transform: returns a TanRaDecToPixel if there are no corrections, or the iterative solver if there are.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
void 
dump(std::ostream &stream) const dumps the transform coefficients to stream.
- 
double 
fit(StarMatchList const &starMatchList) Not implemented yet, because we do it otherwise.
 - 
 
- 
class 
TanRaDecToPixel: public lsst::jointcal::AstrometryTransform - #include <AstrometryTransform.h>
This one is the Tangent Plane (called gnomonic) projection (from celestial sphere to tangent plane)
this transform does not implement corrections, since they are defined the other way around (from pixels to sky), and not invertible analytically. The inversion of tangent point WCS (TanPixelToRaDec) is obtained via inverseTransform().
Public Functions
- 
TanRaDecToPixel(AstrometryTransformLinear const &tan2Pix, Point const &tangentPoint) assume degrees everywhere.
- 
TanRaDecToPixel() 
- 
AstrometryTransformLinear 
getLinPart() const The Linear part (corresponding to CD’s and CRPIX’s)
- 
void 
setTangentPoint(Point const &tangentPoint) Resets the projection (or tangent) point.
- 
Point 
getTangentPoint() const tangent point coordinates (degrees)
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const 
- 
void 
transformPosAndErrors(const FatPoint &in, FatPoint &out) const transform with analytical derivatives
- 
TanPixelToRaDec 
inverted() const exact typed inverse:
- 
std::unique_ptr<AstrometryTransform> 
roughInverse(const Frame ®ion) const Overload the “generic routine” (available for all AstrometryTransform types.
- 
std::unique_ptr<AstrometryTransform> 
inverseTransform(const double precision, const Frame ®ion) const Inverse transform: returns a TanPixelToRaDec.
- 
void 
dump(std::ostream &stream) const dumps the transform coefficients to stream.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
double 
fit(StarMatchList const &starMatchList) fits a transform to a std::list of Point pairs (p1,p2, the Point fields in StarMatch).
After the fit this(p1) yields approximately p2. The returned value is the sum of squared residuals. If you want to fit a partial transform (e.g. such that this(T1(p1)) = T2(p2), use StarMatchList::applyTransform beforehand.
 - 
 
- 
class 
TanSipPixelToRaDec: public lsst::jointcal::BaseTanWcs - #include <AstrometryTransform.h>
Implements the (forward) SIP distorsion scheme.
Public Functions
- 
TanSipPixelToRaDec(AstrometryTransformLinear const &pixToTan, Point const &tangentPoint, const AstrometryTransformPolynomial *corrections = nullptr) pixToTan describes the transform from pix to tangent plane (degrees). TangentPoint in degrees. Corrections are applied before Lin.
- 
AstrometryTransformPolynomial 
getPixelToTangentPlane() const the transformation from pixels to tangent plane (degrees)
- 
virtual void 
pixToTangentPlane(double xPixel, double yPixel, double &xTangentPlane, double &yTangentPlane) const transforms from pixel space to tangent plane (degrees)
- 
TanSipPixelToRaDec() 
- 
std::unique_ptr<AstrometryTransform> 
inverseTransform(const double precision, const Frame ®ion) const Inverse transform: returns a TanRaDecToPixel if there are no corrections, or the iterative solver if there are.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
- 
void 
dump(std::ostream &stream) const dumps the transform coefficients to stream.
- 
double 
fit(StarMatchList const &starMatchList) Not implemented yet, because we do it otherwise.
 - 
 
- 
class 
UserTransform: public lsst::jointcal::AstrometryTransform - #include <AstrometryTransform.h>
A run-time transform that allows users to define a AstrometryTransform with minimal coding (just the apply routine).
Public Functions
- 
UserTransform(AstrometryTransformFun &fun, const void *userData) the transform routine and extra data that it may need.
- 
void 
apply(const double xIn, const double yIn, double &xOut, double &yOut) const 
- 
double 
fit(StarMatchList const &starMatchList) fits a transform to a std::list of Point pairs (p1,p2, the Point fields in StarMatch).
After the fit this(p1) yields approximately p2. The returned value is the sum of squared residuals. If you want to fit a partial transform (e.g. such that this(T1(p1)) = T2(p2), use StarMatchList::applyTransform beforehand.
- 
std::unique_ptr<AstrometryTransform> 
clone() const returns a copy (allocated by new) of the transformation.
 - 
 
- 
typedef